1be1d678aSKris Buschelman 2b97cf49bSBarry Smith /* 3b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 4b97cf49bSBarry Smith */ 5e37d9405SJed Brown #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/ 6841d17a1SFande Kong #include <petscsf.h> 7b97cf49bSBarry Smith 840244768SBarry Smith /* 97dae84e0SHong Zhang * The interface should be easy to use for both MatCreateSubMatrix (parallel sub-matrix) and MatCreateSubMatrices (sequential sub-matrices) 1040244768SBarry Smith * */ 117dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_MPIAdj_data(Mat adj,IS irows, IS icols, PetscInt **sadj_xadj,PetscInt **sadj_adjncy,PetscInt **sadj_values) 1240244768SBarry Smith { 13131c27b5Sprj- PetscInt nlrows_is,icols_n,i,j,nroots,nleaves,rlocalindex,*ncols_send,*ncols_recv; 1440244768SBarry Smith PetscInt nlrows_mat,*adjncy_recv,Ncols_recv,Ncols_send,*xadj_recv,*values_recv; 15e895ccc0SFande Kong PetscInt *ncols_recv_offsets,loc,rnclos,*sadjncy,*sxadj,*svalues; 1640244768SBarry Smith const PetscInt *irows_indices,*icols_indices,*xadj, *adjncy; 17131c27b5Sprj- PetscMPIInt owner; 1840244768SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)adj->data; 1940244768SBarry Smith PetscLayout rmap; 2040244768SBarry Smith MPI_Comm comm; 2140244768SBarry Smith PetscSF sf; 2240244768SBarry Smith PetscSFNode *iremote; 2340244768SBarry Smith PetscBool done; 2440244768SBarry Smith 2540244768SBarry Smith PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)adj,&comm)); 279566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(adj,&rmap,NULL)); 289566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irows,&nlrows_is)); 299566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irows,&irows_indices)); 309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlrows_is,&iremote)); 3140244768SBarry Smith /* construct sf graph*/ 3240244768SBarry Smith nleaves = nlrows_is; 3340244768SBarry Smith for (i=0; i<nlrows_is; i++) { 3440244768SBarry Smith owner = -1; 3540244768SBarry Smith rlocalindex = -1; 369566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwnerIndex(rmap,irows_indices[i],&owner,&rlocalindex)); 3740244768SBarry Smith iremote[i].rank = owner; 3840244768SBarry Smith iremote[i].index = rlocalindex; 3940244768SBarry Smith } 409566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done)); 419566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(nlrows_mat,&ncols_send,nlrows_is,&xadj_recv,nlrows_is+1,&ncols_recv_offsets,nlrows_is,&ncols_recv)); 4240244768SBarry Smith nroots = nlrows_mat; 4340244768SBarry Smith for (i=0; i<nlrows_mat; i++) { 4440244768SBarry Smith ncols_send[i] = xadj[i+1]-xadj[i]; 4540244768SBarry Smith } 469566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm,&sf)); 479566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER)); 489566063dSJacob Faibussowitsch PetscCall(PetscSFSetType(sf,PETSCSFBASIC)); 499566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf)); 509566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv,MPI_REPLACE)); 519566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv,MPI_REPLACE)); 529566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv,MPI_REPLACE)); 539566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv,MPI_REPLACE)); 549566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 5540244768SBarry Smith Ncols_recv =0; 5640244768SBarry Smith for (i=0; i<nlrows_is; i++) { 5740244768SBarry Smith Ncols_recv += ncols_recv[i]; 5840244768SBarry Smith ncols_recv_offsets[i+1] = ncols_recv[i]+ncols_recv_offsets[i]; 5940244768SBarry Smith } 6040244768SBarry Smith Ncols_send = 0; 6140244768SBarry Smith for (i=0; i<nlrows_mat; i++) { 6240244768SBarry Smith Ncols_send += ncols_send[i]; 6340244768SBarry Smith } 649566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Ncols_recv,&iremote)); 659566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Ncols_recv,&adjncy_recv)); 6640244768SBarry Smith nleaves = Ncols_recv; 6740244768SBarry Smith Ncols_recv = 0; 6840244768SBarry Smith for (i=0; i<nlrows_is; i++) { 699566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(rmap,irows_indices[i],&owner)); 7040244768SBarry Smith for (j=0; j<ncols_recv[i]; j++) { 7140244768SBarry Smith iremote[Ncols_recv].rank = owner; 7240244768SBarry Smith iremote[Ncols_recv++].index = xadj_recv[i]+j; 7340244768SBarry Smith } 7440244768SBarry Smith } 759566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irows,&irows_indices)); 7640244768SBarry Smith /*if we need to deal with edge weights ???*/ 779566063dSJacob Faibussowitsch if (a->useedgeweights) PetscCall(PetscCalloc1(Ncols_recv,&values_recv)); 7840244768SBarry Smith nroots = Ncols_send; 799566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm,&sf)); 809566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER)); 819566063dSJacob Faibussowitsch PetscCall(PetscSFSetType(sf,PETSCSFBASIC)); 829566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf)); 839566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv,MPI_REPLACE)); 849566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv,MPI_REPLACE)); 85e895ccc0SFande Kong if (a->useedgeweights) { 869566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv,MPI_REPLACE)); 879566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv,MPI_REPLACE)); 8840244768SBarry Smith } 899566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 909566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done)); 919566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icols,&icols_n)); 929566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icols,&icols_indices)); 9340244768SBarry Smith rnclos = 0; 9440244768SBarry Smith for (i=0; i<nlrows_is; i++) { 9540244768SBarry Smith for (j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++) { 969566063dSJacob Faibussowitsch PetscCall(PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc)); 9740244768SBarry Smith if (loc<0) { 9840244768SBarry Smith adjncy_recv[j] = -1; 99e895ccc0SFande Kong if (a->useedgeweights) values_recv[j] = -1; 10040244768SBarry Smith ncols_recv[i]--; 10140244768SBarry Smith } else { 10240244768SBarry Smith rnclos++; 10340244768SBarry Smith } 10440244768SBarry Smith } 10540244768SBarry Smith } 1069566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icols,&icols_indices)); 1079566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(rnclos,&sadjncy)); 1089566063dSJacob Faibussowitsch if (a->useedgeweights) PetscCall(PetscCalloc1(rnclos,&svalues)); 1099566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nlrows_is+1,&sxadj)); 11040244768SBarry Smith rnclos = 0; 11140244768SBarry Smith for (i=0; i<nlrows_is; i++) { 11240244768SBarry Smith for (j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++) { 11340244768SBarry Smith if (adjncy_recv[j]<0) continue; 11440244768SBarry Smith sadjncy[rnclos] = adjncy_recv[j]; 115e895ccc0SFande Kong if (a->useedgeweights) svalues[rnclos] = values_recv[j]; 11640244768SBarry Smith rnclos++; 11740244768SBarry Smith } 11840244768SBarry Smith } 11940244768SBarry Smith for (i=0; i<nlrows_is; i++) { 12040244768SBarry Smith sxadj[i+1] = sxadj[i]+ncols_recv[i]; 12140244768SBarry Smith } 1229566063dSJacob Faibussowitsch if (sadj_xadj) { *sadj_xadj = sxadj;} else PetscCall(PetscFree(sxadj)); 1239566063dSJacob Faibussowitsch if (sadj_adjncy) { *sadj_adjncy = sadjncy;} else PetscCall(PetscFree(sadjncy)); 12440244768SBarry Smith if (sadj_values) { 125f4259b30SLisandro Dalcin if (a->useedgeweights) *sadj_values = svalues; else *sadj_values=NULL; 12640244768SBarry Smith } else { 1279566063dSJacob Faibussowitsch if (a->useedgeweights) PetscCall(PetscFree(svalues)); 12840244768SBarry Smith } 1299566063dSJacob Faibussowitsch PetscCall(PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(adjncy_recv)); 1319566063dSJacob Faibussowitsch if (a->useedgeweights) PetscCall(PetscFree(values_recv)); 13240244768SBarry Smith PetscFunctionReturn(0); 13340244768SBarry Smith } 13440244768SBarry Smith 1357dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj_Private(Mat mat,PetscInt n,const IS irow[],const IS icol[],PetscBool subcomm,MatReuse scall,Mat *submat[]) 13640244768SBarry Smith { 13740244768SBarry Smith PetscInt i,irow_n,icol_n,*sxadj,*sadjncy,*svalues; 13840244768SBarry Smith PetscInt *indices,nindx,j,k,loc; 13940244768SBarry Smith PetscMPIInt issame; 14040244768SBarry Smith const PetscInt *irow_indices,*icol_indices; 14140244768SBarry Smith MPI_Comm scomm_row,scomm_col,scomm_mat; 14240244768SBarry Smith 14340244768SBarry Smith PetscFunctionBegin; 14440244768SBarry Smith nindx = 0; 14540244768SBarry Smith /* 14640244768SBarry Smith * Estimate a maximum number for allocating memory 14740244768SBarry Smith */ 14840244768SBarry Smith for (i=0; i<n; i++) { 1499566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i],&irow_n)); 1509566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i],&icol_n)); 15140244768SBarry Smith nindx = nindx>(irow_n+icol_n)? nindx:(irow_n+icol_n); 15240244768SBarry Smith } 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nindx,&indices)); 15440244768SBarry Smith /* construct a submat */ 15540244768SBarry Smith for (i=0; i<n; i++) { 15640244768SBarry Smith if (subcomm) { 1579566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)irow[i],&scomm_row)); 1589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)icol[i],&scomm_col)); 1599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row,scomm_col,&issame)); 16008401ef6SPierre Jolivet PetscCheck(issame == MPI_IDENT,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row index set must have the same comm as the col index set"); 1619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame)); 16208401ef6SPierre Jolivet PetscCheck(issame != MPI_IDENT,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP," can not use PETSC_COMM_SELF as comm when extracting a parallel submatrix"); 16340244768SBarry Smith } else { 16440244768SBarry Smith scomm_row = PETSC_COMM_SELF; 16540244768SBarry Smith } 16640244768SBarry Smith /*get sub-matrix data*/ 167f4259b30SLisandro Dalcin sxadj=NULL; sadjncy=NULL; svalues=NULL; 1689566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues)); 1699566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i],&irow_n)); 1709566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i],&icol_n)); 1719566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irow[i],&irow_indices)); 1729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indices,irow_indices,irow_n)); 1739566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i],&irow_indices)); 1749566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icol[i],&icol_indices)); 1759566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indices+irow_n,icol_indices,icol_n)); 1769566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i],&icol_indices)); 17740244768SBarry Smith nindx = irow_n+icol_n; 1789566063dSJacob Faibussowitsch PetscCall(PetscSortRemoveDupsInt(&nindx,indices)); 17940244768SBarry Smith /* renumber columns */ 18040244768SBarry Smith for (j=0; j<irow_n; j++) { 18140244768SBarry Smith for (k=sxadj[j]; k<sxadj[j+1]; k++) { 1829566063dSJacob Faibussowitsch PetscCall(PetscFindInt(sadjncy[k],nindx,indices,&loc)); 18308401ef6SPierre Jolivet PetscCheck(loc>=0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %" PetscInt_FMT,sadjncy[k]); 18440244768SBarry Smith sadjncy[k] = loc; 18540244768SBarry Smith } 18640244768SBarry Smith } 18740244768SBarry Smith if (scall==MAT_INITIAL_MATRIX) { 1889566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i])); 18940244768SBarry Smith } else { 19040244768SBarry Smith Mat sadj = *(submat[i]); 19140244768SBarry Smith Mat_MPIAdj *sa = (Mat_MPIAdj*)((sadj)->data); 1929566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sadj,&scomm_mat)); 1939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row,scomm_mat,&issame)); 19408401ef6SPierre Jolivet PetscCheck(issame == MPI_IDENT,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"submatrix must have the same comm as the col index set"); 1959566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sa->i,sxadj,irow_n+1)); 1969566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sa->j,sadjncy,sxadj[irow_n])); 1979566063dSJacob Faibussowitsch if (svalues) PetscCall(PetscArraycpy(sa->values,svalues,sxadj[irow_n])); 1989566063dSJacob Faibussowitsch PetscCall(PetscFree(sxadj)); 1999566063dSJacob Faibussowitsch PetscCall(PetscFree(sadjncy)); 2009566063dSJacob Faibussowitsch if (svalues) PetscCall(PetscFree(svalues)); 20140244768SBarry Smith } 20240244768SBarry Smith } 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 20440244768SBarry Smith PetscFunctionReturn(0); 20540244768SBarry Smith } 20640244768SBarry Smith 2077dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatricesMPI_MPIAdj(Mat mat,PetscInt n, const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 20840244768SBarry Smith { 20940244768SBarry Smith /*get sub-matrices across a sub communicator */ 21040244768SBarry Smith PetscFunctionBegin; 2119566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat)); 21240244768SBarry Smith PetscFunctionReturn(0); 21340244768SBarry Smith } 21440244768SBarry Smith 2157dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 21640244768SBarry Smith { 21740244768SBarry Smith PetscFunctionBegin; 21840244768SBarry Smith /*get sub-matrices based on PETSC_COMM_SELF */ 2199566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat)); 22040244768SBarry Smith PetscFunctionReturn(0); 22140244768SBarry Smith } 22240244768SBarry Smith 22340244768SBarry Smith static PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 224b97cf49bSBarry Smith { 2253eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 226d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 2272dcb1b2aSMatthew Knepley const char *name; 228ce1411ecSBarry Smith PetscViewerFormat format; 229b97cf49bSBarry Smith 230433994e6SBarry Smith PetscFunctionBegin; 2319566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)A,&name)); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 233fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 2343a40ed3dSBarry Smith PetscFunctionReturn(0); 235f7d195e4SLawrence Mitchell } else { 236f7d195e4SLawrence Mitchell PetscCheck(format != PETSC_VIEWER_ASCII_MATLAB,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported"); 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 239b97cf49bSBarry Smith for (i=0; i<m; i++) { 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"row %" PetscInt_FMT ":",i+A->rmap->rstart)); 241b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 2420170919eSFande Kong if (a->values) { 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," (%" PetscInt_FMT ", %" PetscInt_FMT ") ",a->j[j], a->values[j])); 2440170919eSFande Kong } else { 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT " ",a->j[j])); 2460752156aSBarry Smith } 2470170919eSFande Kong } 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n")); 249b97cf49bSBarry Smith } 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); 2519566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2537b23a99aSBarry Smith } 2543a40ed3dSBarry Smith PetscFunctionReturn(0); 255b97cf49bSBarry Smith } 256b97cf49bSBarry Smith 25740244768SBarry Smith static PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 258b97cf49bSBarry Smith { 259ace3abfcSBarry Smith PetscBool iascii; 260b97cf49bSBarry Smith 261433994e6SBarry Smith PetscFunctionBegin; 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 2631baa6e33SBarry Smith if (iascii) PetscCall(MatView_MPIAdj_ASCII(A,viewer)); 2643a40ed3dSBarry Smith PetscFunctionReturn(0); 265b97cf49bSBarry Smith } 266b97cf49bSBarry Smith 26740244768SBarry Smith static PetscErrorCode MatDestroy_MPIAdj(Mat mat) 268b97cf49bSBarry Smith { 2693eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 27094d884c6SBarry Smith 27194d884c6SBarry Smith PetscFunctionBegin; 272aa482453SBarry Smith #if defined(PETSC_USE_LOG) 273c0aa6a63SJacob Faibussowitsch PetscLogObjectState((PetscObject)mat,"Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT,mat->rmap->n,mat->cmap->n,a->nz); 274b97cf49bSBarry Smith #endif 2759566063dSJacob Faibussowitsch PetscCall(PetscFree(a->diag)); 2760400557aSBarry Smith if (a->freeaij) { 27714196391SBarry Smith if (a->freeaijwithfree) { 27814196391SBarry Smith if (a->i) free(a->i); 27914196391SBarry Smith if (a->j) free(a->j); 28014196391SBarry Smith } else { 2819566063dSJacob Faibussowitsch PetscCall(PetscFree(a->i)); 2829566063dSJacob Faibussowitsch PetscCall(PetscFree(a->j)); 2839566063dSJacob Faibussowitsch PetscCall(PetscFree(a->values)); 28414196391SBarry Smith } 2850400557aSBarry Smith } 2869566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 2879566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 2889566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)mat,NULL)); 2899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL)); 2909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL)); 2912e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjToSeq_C",NULL)); 2923a40ed3dSBarry Smith PetscFunctionReturn(0); 293b97cf49bSBarry Smith } 294b97cf49bSBarry Smith 29540244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 296b97cf49bSBarry Smith { 2973eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 298b97cf49bSBarry Smith 299433994e6SBarry Smith PetscFunctionBegin; 30012c028f9SKris Buschelman switch (op) { 30177e54ba9SKris Buschelman case MAT_SYMMETRIC: 30212c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3039a4540c5SBarry Smith case MAT_HERMITIAN: 3044e0d8c25SBarry Smith a->symmetric = flg; 3059a4540c5SBarry Smith break; 3069a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 3079a4540c5SBarry Smith break; 30812c028f9SKris Buschelman default: 3099566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"Option %s ignored\n",MatOptions[op])); 31012c028f9SKris Buschelman break; 311b97cf49bSBarry Smith } 3123a40ed3dSBarry Smith PetscFunctionReturn(0); 313b97cf49bSBarry Smith } 314b97cf49bSBarry Smith 31540244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 316b97cf49bSBarry Smith { 3173eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 318b97cf49bSBarry Smith 319433994e6SBarry Smith PetscFunctionBegin; 320d0f46423SBarry Smith row -= A->rmap->rstart; 321aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 322b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 3232b1d8763SJed Brown if (v) { 3242b1d8763SJed Brown PetscInt j; 3252b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 3269566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 3272b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz); 3289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->rowvalues_alloc,&a->rowvalues)); 329b97cf49bSBarry Smith } 33091f8cafeSFande Kong for (j=0; j<*nz; j++) { 33191f8cafeSFande Kong a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0; 33291f8cafeSFande Kong } 3332b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 334b97cf49bSBarry Smith } 33592b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 3363a40ed3dSBarry Smith PetscFunctionReturn(0); 337b97cf49bSBarry Smith } 338b97cf49bSBarry Smith 33940244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 340b97cf49bSBarry Smith { 341433994e6SBarry Smith PetscFunctionBegin; 3423a40ed3dSBarry Smith PetscFunctionReturn(0); 343b97cf49bSBarry Smith } 344b97cf49bSBarry Smith 34540244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 346b97cf49bSBarry Smith { 3473eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 348ace3abfcSBarry Smith PetscBool flag; 349b97cf49bSBarry Smith 350433994e6SBarry Smith PetscFunctionBegin; 351b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 352d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 3530f5bd95cSBarry Smith flag = PETSC_FALSE; 354b97cf49bSBarry Smith } 355b97cf49bSBarry Smith 356b97cf49bSBarry Smith /* if the a->i are the same */ 3579566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(a->i,b->i,A->rmap->n+1,&flag)); 358b97cf49bSBarry Smith 359b97cf49bSBarry Smith /* if a->j are the same */ 3609566063dSJacob Faibussowitsch PetscCall(PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag)); 361b97cf49bSBarry Smith 3621c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 3633a40ed3dSBarry Smith PetscFunctionReturn(0); 364b97cf49bSBarry Smith } 365b97cf49bSBarry Smith 36640244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 3676cd91f64SBarry Smith { 368b24ad042SBarry Smith PetscInt i; 3696cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3701a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 3716cd91f64SBarry Smith 3726cd91f64SBarry Smith PetscFunctionBegin; 373d0f46423SBarry Smith *m = A->rmap->n; 3746cd91f64SBarry Smith *ia = a->i; 3756cd91f64SBarry Smith *ja = a->j; 3766cd91f64SBarry Smith *done = PETSC_TRUE; 377d42735a1SBarry Smith if (oshift) { 378d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 379d42735a1SBarry Smith (*ja)[i]++; 380d42735a1SBarry Smith } 381d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 382d42735a1SBarry Smith } 383d42735a1SBarry Smith PetscFunctionReturn(0); 384d42735a1SBarry Smith } 385d42735a1SBarry Smith 38640244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 387d42735a1SBarry Smith { 388b24ad042SBarry Smith PetscInt i; 389d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3901a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 391d42735a1SBarry Smith 392d42735a1SBarry Smith PetscFunctionBegin; 393aed4548fSBarry Smith PetscCheck(!ia || a->i == *ia,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 394aed4548fSBarry Smith PetscCheck(!ja || a->j == *ja,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 395d42735a1SBarry Smith if (oshift) { 39628b400f6SJacob Faibussowitsch PetscCheck(ia,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument"); 39728b400f6SJacob Faibussowitsch PetscCheck(ja,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument"); 398d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 399d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 400d42735a1SBarry Smith (*ja)[i]--; 401d42735a1SBarry Smith } 402d42735a1SBarry Smith } 4036cd91f64SBarry Smith PetscFunctionReturn(0); 4046cd91f64SBarry Smith } 405b97cf49bSBarry Smith 40619fd82e9SBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 40717667f90SBarry Smith { 40817667f90SBarry Smith Mat B; 40917667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 41017667f90SBarry Smith const PetscInt *rj; 41117667f90SBarry Smith const PetscScalar *ra; 41217667f90SBarry Smith MPI_Comm comm; 41317667f90SBarry Smith 41417667f90SBarry Smith PetscFunctionBegin; 4159566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 4169566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 4179566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,NULL)); 41817667f90SBarry Smith 41917667f90SBarry Smith /* count the number of nonzeros per row */ 42017667f90SBarry Smith for (i=0; i<m; i++) { 4219566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i+rstart,&len,&rj,NULL)); 4225ee9ba1cSJed Brown for (j=0; j<len; j++) { 4235ee9ba1cSJed Brown if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 4245ee9ba1cSJed Brown } 42517667f90SBarry Smith nzeros += len; 4269566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i+rstart,&len,&rj,NULL)); 42717667f90SBarry Smith } 42817667f90SBarry Smith 42917667f90SBarry Smith /* malloc space for nonzeros */ 4309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros+1,&a)); 4319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N+1,&ia)); 4329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros+1,&ja)); 43317667f90SBarry Smith 43417667f90SBarry Smith nzeros = 0; 43517667f90SBarry Smith ia[0] = 0; 43617667f90SBarry Smith for (i=0; i<m; i++) { 4379566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i+rstart,&len,&rj,&ra)); 43817667f90SBarry Smith cnt = 0; 43917667f90SBarry Smith for (j=0; j<len; j++) { 4405ee9ba1cSJed Brown if (rj[j] != i+rstart) { /* if not diagonal */ 44117667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 44217667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 44317667f90SBarry Smith } 4445ee9ba1cSJed Brown } 4459566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i+rstart,&len,&rj,&ra)); 44617667f90SBarry Smith nzeros += cnt; 44717667f90SBarry Smith ia[i+1] = nzeros; 44817667f90SBarry Smith } 44917667f90SBarry Smith 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A,&comm)); 4519566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&B)); 4529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N)); 4539566063dSJacob Faibussowitsch PetscCall(MatSetType(B,type)); 4549566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(B,ia,ja,a)); 45517667f90SBarry Smith 456511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 4579566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A,&B)); 45817667f90SBarry Smith } else { 45917667f90SBarry Smith *newmat = B; 46017667f90SBarry Smith } 46117667f90SBarry Smith PetscFunctionReturn(0); 46217667f90SBarry Smith } 46317667f90SBarry Smith 464b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 465f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL, 4663eda8832SBarry Smith MatGetRow_MPIAdj, 4673eda8832SBarry Smith MatRestoreRow_MPIAdj, 468f4259b30SLisandro Dalcin NULL, 469f4259b30SLisandro Dalcin /* 4*/ NULL, 470f4259b30SLisandro Dalcin NULL, 471f4259b30SLisandro Dalcin NULL, 472f4259b30SLisandro Dalcin NULL, 473f4259b30SLisandro Dalcin NULL, 474f4259b30SLisandro Dalcin NULL, 475f4259b30SLisandro Dalcin /*10*/ NULL, 476f4259b30SLisandro Dalcin NULL, 477f4259b30SLisandro Dalcin NULL, 478f4259b30SLisandro Dalcin NULL, 479f4259b30SLisandro Dalcin NULL, 480f4259b30SLisandro Dalcin /*15*/ NULL, 4813eda8832SBarry Smith MatEqual_MPIAdj, 482f4259b30SLisandro Dalcin NULL, 483f4259b30SLisandro Dalcin NULL, 484f4259b30SLisandro Dalcin NULL, 485f4259b30SLisandro Dalcin /*20*/ NULL, 486f4259b30SLisandro Dalcin NULL, 4873eda8832SBarry Smith MatSetOption_MPIAdj, 488f4259b30SLisandro Dalcin NULL, 489f4259b30SLisandro Dalcin /*24*/ NULL, 490f4259b30SLisandro Dalcin NULL, 491f4259b30SLisandro Dalcin NULL, 492f4259b30SLisandro Dalcin NULL, 493f4259b30SLisandro Dalcin NULL, 494f4259b30SLisandro Dalcin /*29*/ NULL, 495f4259b30SLisandro Dalcin NULL, 496f4259b30SLisandro Dalcin NULL, 497f4259b30SLisandro Dalcin NULL, 498f4259b30SLisandro Dalcin NULL, 499f4259b30SLisandro Dalcin /*34*/ NULL, 500f4259b30SLisandro Dalcin NULL, 501f4259b30SLisandro Dalcin NULL, 502f4259b30SLisandro Dalcin NULL, 503f4259b30SLisandro Dalcin NULL, 504f4259b30SLisandro Dalcin /*39*/ NULL, 5057dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 506f4259b30SLisandro Dalcin NULL, 507f4259b30SLisandro Dalcin NULL, 508f4259b30SLisandro Dalcin NULL, 509f4259b30SLisandro Dalcin /*44*/ NULL, 510f4259b30SLisandro Dalcin NULL, 5117d68702bSBarry Smith MatShift_Basic, 512f4259b30SLisandro Dalcin NULL, 513f4259b30SLisandro Dalcin NULL, 514f4259b30SLisandro Dalcin /*49*/ NULL, 5156cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 516d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 517f4259b30SLisandro Dalcin NULL, 518f4259b30SLisandro Dalcin NULL, 519f4259b30SLisandro Dalcin /*54*/ NULL, 520f4259b30SLisandro Dalcin NULL, 521f4259b30SLisandro Dalcin NULL, 522f4259b30SLisandro Dalcin NULL, 523f4259b30SLisandro Dalcin NULL, 524f4259b30SLisandro Dalcin /*59*/ NULL, 525b9b97703SBarry Smith MatDestroy_MPIAdj, 526b9b97703SBarry Smith MatView_MPIAdj, 52717667f90SBarry Smith MatConvertFrom_MPIAdj, 528f4259b30SLisandro Dalcin NULL, 529f4259b30SLisandro Dalcin /*64*/ NULL, 530f4259b30SLisandro Dalcin NULL, 531f4259b30SLisandro Dalcin NULL, 532f4259b30SLisandro Dalcin NULL, 533f4259b30SLisandro Dalcin NULL, 534f4259b30SLisandro Dalcin /*69*/ NULL, 535f4259b30SLisandro Dalcin NULL, 536f4259b30SLisandro Dalcin NULL, 537f4259b30SLisandro Dalcin NULL, 538f4259b30SLisandro Dalcin NULL, 539f4259b30SLisandro Dalcin /*74*/ NULL, 540f4259b30SLisandro Dalcin NULL, 541f4259b30SLisandro Dalcin NULL, 542f4259b30SLisandro Dalcin NULL, 543f4259b30SLisandro Dalcin NULL, 544f4259b30SLisandro Dalcin /*79*/ NULL, 545f4259b30SLisandro Dalcin NULL, 546f4259b30SLisandro Dalcin NULL, 547f4259b30SLisandro Dalcin NULL, 548f4259b30SLisandro Dalcin NULL, 549f4259b30SLisandro Dalcin /*84*/ NULL, 550f4259b30SLisandro Dalcin NULL, 551f4259b30SLisandro Dalcin NULL, 552f4259b30SLisandro Dalcin NULL, 553f4259b30SLisandro Dalcin NULL, 554f4259b30SLisandro Dalcin /*89*/ NULL, 555f4259b30SLisandro Dalcin NULL, 556f4259b30SLisandro Dalcin NULL, 557f4259b30SLisandro Dalcin NULL, 558f4259b30SLisandro Dalcin NULL, 559f4259b30SLisandro Dalcin /*94*/ NULL, 560f4259b30SLisandro Dalcin NULL, 561f4259b30SLisandro Dalcin NULL, 562f4259b30SLisandro Dalcin NULL, 563f4259b30SLisandro Dalcin NULL, 564f4259b30SLisandro Dalcin /*99*/ NULL, 565f4259b30SLisandro Dalcin NULL, 566f4259b30SLisandro Dalcin NULL, 567f4259b30SLisandro Dalcin NULL, 568f4259b30SLisandro Dalcin NULL, 569f4259b30SLisandro Dalcin /*104*/ NULL, 570f4259b30SLisandro Dalcin NULL, 571f4259b30SLisandro Dalcin NULL, 572f4259b30SLisandro Dalcin NULL, 573f4259b30SLisandro Dalcin NULL, 574f4259b30SLisandro Dalcin /*109*/ NULL, 575f4259b30SLisandro Dalcin NULL, 576f4259b30SLisandro Dalcin NULL, 577f4259b30SLisandro Dalcin NULL, 578f4259b30SLisandro Dalcin NULL, 579f4259b30SLisandro Dalcin /*114*/ NULL, 580f4259b30SLisandro Dalcin NULL, 581f4259b30SLisandro Dalcin NULL, 582f4259b30SLisandro Dalcin NULL, 583f4259b30SLisandro Dalcin NULL, 584f4259b30SLisandro Dalcin /*119*/ NULL, 585f4259b30SLisandro Dalcin NULL, 586f4259b30SLisandro Dalcin NULL, 587f4259b30SLisandro Dalcin NULL, 588f4259b30SLisandro Dalcin NULL, 589f4259b30SLisandro Dalcin /*124*/ NULL, 590f4259b30SLisandro Dalcin NULL, 591f4259b30SLisandro Dalcin NULL, 592f4259b30SLisandro Dalcin NULL, 5937dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 594f4259b30SLisandro Dalcin /*129*/ NULL, 595f4259b30SLisandro Dalcin NULL, 596f4259b30SLisandro Dalcin NULL, 597f4259b30SLisandro Dalcin NULL, 598f4259b30SLisandro Dalcin NULL, 599f4259b30SLisandro Dalcin /*134*/ NULL, 600f4259b30SLisandro Dalcin NULL, 601f4259b30SLisandro Dalcin NULL, 602f4259b30SLisandro Dalcin NULL, 603f4259b30SLisandro Dalcin NULL, 604f4259b30SLisandro Dalcin /*139*/ NULL, 605f4259b30SLisandro Dalcin NULL, 606d70f29a3SPierre Jolivet NULL, 607d70f29a3SPierre Jolivet NULL, 608d70f29a3SPierre Jolivet NULL, 609d70f29a3SPierre Jolivet /*144*/NULL, 610d70f29a3SPierre Jolivet NULL, 611d70f29a3SPierre Jolivet NULL, 612*99a7f59eSMark Adams NULL, 613*99a7f59eSMark Adams NULL, 614f4259b30SLisandro Dalcin NULL 6153964eb88SJed Brown }; 616b97cf49bSBarry Smith 617f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 618a23d5eceSKris Buschelman { 619a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 620e895ccc0SFande Kong PetscBool useedgeweights; 621a23d5eceSKris Buschelman 622a23d5eceSKris Buschelman PetscFunctionBegin; 6239566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->rmap)); 6249566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->cmap)); 625e895ccc0SFande Kong if (values) useedgeweights = PETSC_TRUE; else useedgeweights = PETSC_FALSE; 626e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */ 6271c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((int*)&useedgeweights,(int*)&b->useedgeweights,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)B))); 62858ec18a5SLisandro Dalcin 62976bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 63076bd3646SJed Brown PetscInt ii; 63176bd3646SJed Brown 63208401ef6SPierre Jolivet PetscCheck(i[0] == 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %" PetscInt_FMT,i[0]); 633d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 634aed4548fSBarry Smith PetscCheck(i[ii] >= 0 && i[ii] >= i[ii-1],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i[%" PetscInt_FMT "]=%" PetscInt_FMT " index is out of range: i[%" PetscInt_FMT "]=%" PetscInt_FMT,ii,i[ii],ii-1,i[ii-1]); 635a23d5eceSKris Buschelman } 636d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 637aed4548fSBarry Smith PetscCheck(j[ii] >= 0 && j[ii] < B->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %" PetscInt_FMT " out of range %" PetscInt_FMT,ii,j[ii]); 638a23d5eceSKris Buschelman } 63976bd3646SJed Brown } 64058ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 641a23d5eceSKris Buschelman 642a23d5eceSKris Buschelman b->j = j; 643a23d5eceSKris Buschelman b->i = i; 644a23d5eceSKris Buschelman b->values = values; 645a23d5eceSKris Buschelman 646d0f46423SBarry Smith b->nz = i[B->rmap->n]; 647f4259b30SLisandro Dalcin b->diag = NULL; 648a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 649a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 650a23d5eceSKris Buschelman 6519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 6529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 653a23d5eceSKris Buschelman PetscFunctionReturn(0); 654a23d5eceSKris Buschelman } 655b97cf49bSBarry Smith 656f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 6579a3665c6SJed Brown { 6589a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 6599a3665c6SJed Brown const PetscInt *ranges; 6609a3665c6SJed Brown MPI_Comm acomm,bcomm; 6619a3665c6SJed Brown MPI_Group agroup,bgroup; 6629a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 6639a3665c6SJed Brown 6649a3665c6SJed Brown PetscFunctionBegin; 6650298fd71SBarry Smith *B = NULL; 6669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A,&acomm)); 6679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm,&size)); 6689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm,&rank)); 6699566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A,&ranges)); 6709a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6719a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 6729a3665c6SJed Brown } 6739a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 6749566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 6759a3665c6SJed Brown *B = A; 6769a3665c6SJed Brown PetscFunctionReturn(0); 6779a3665c6SJed Brown } 6789a3665c6SJed Brown 6799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nranks,&ranks)); 6809a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6819a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 6829a3665c6SJed Brown } 6839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(acomm,&agroup)); 6849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(agroup,nranks,ranks,&bgroup)); 6859566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks)); 6869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(acomm,bgroup,&bcomm)); 6879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&agroup)); 6889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&bgroup)); 6899a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 6909a3665c6SJed Brown PetscInt m,N; 6919a3665c6SJed Brown Mat_MPIAdj *b; 6929566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 6939566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 6949566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B)); 6959a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 6969a3665c6SJed Brown b->freeaij = PETSC_FALSE; 6979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&bcomm)); 6989a3665c6SJed Brown } 6999a3665c6SJed Brown PetscFunctionReturn(0); 7009a3665c6SJed Brown } 7019a3665c6SJed Brown 702b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B) 703b44e7856SBarry Smith { 704b44e7856SBarry Smith PetscInt M,N,*II,*J,NZ,nz,m,nzstart,i; 705b44e7856SBarry Smith PetscInt *Values = NULL; 706b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 707b44e7856SBarry Smith PetscMPIInt mnz,mm,*allnz,*allm,size,*dispnz,*dispm; 708b44e7856SBarry Smith 709b44e7856SBarry Smith PetscFunctionBegin; 7109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 7119566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 7129566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 713b44e7856SBarry Smith nz = adj->nz; 71408401ef6SPierre Jolivet PetscCheck(adj->i[m] == nz,PETSC_COMM_SELF,PETSC_ERR_PLIB,"nz %" PetscInt_FMT " not correct i[m] %" PetscInt_FMT,nz,adj->i[m]); 7159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A))); 716d4e69552SBarry Smith 7179566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nz,&mnz)); 7189566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size,&allnz,size,&dispnz)); 7199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A))); 720b44e7856SBarry Smith dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1]; 721b44e7856SBarry Smith if (adj->values) { 7229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ,&Values)); 7239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A))); 724b44e7856SBarry Smith } 7259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ,&J)); 7269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A))); 7279566063dSJacob Faibussowitsch PetscCall(PetscFree2(allnz,dispnz)); 7289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A))); 729b44e7856SBarry Smith nzstart -= nz; 730b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 731b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] += nzstart; 7329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M+1,&II)); 7339566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(m,&mm)); 7349566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size,&allm,size,&dispm)); 7359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A))); 736b44e7856SBarry Smith dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1]; 7379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A))); 7389566063dSJacob Faibussowitsch PetscCall(PetscFree2(allm,dispm)); 739b44e7856SBarry Smith II[M] = NZ; 740b44e7856SBarry Smith /* shift the i[] values back */ 741b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] -= nzstart; 7429566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B)); 743b44e7856SBarry Smith PetscFunctionReturn(0); 744b44e7856SBarry Smith } 745b44e7856SBarry Smith 7469a3665c6SJed Brown /*@ 7479a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 7489a3665c6SJed Brown 7499a3665c6SJed Brown Collective 7509a3665c6SJed Brown 7514165533cSJose E. Roman Input Parameter: 7529a3665c6SJed Brown . A - original MPIAdj matrix 7539a3665c6SJed Brown 7544165533cSJose E. Roman Output Parameter: 7550298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 7569a3665c6SJed Brown 7579a3665c6SJed Brown Level: developer 7589a3665c6SJed Brown 7599a3665c6SJed Brown Note: 7609a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 7619a3665c6SJed Brown 7629a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 7639a3665c6SJed Brown 764db781477SPatrick Sanan .seealso: `MatCreateMPIAdj()` 7659a3665c6SJed Brown @*/ 7669a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 7679a3665c6SJed Brown { 7689a3665c6SJed Brown PetscFunctionBegin; 7699a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 770cac4c232SBarry Smith PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B)); 7719a3665c6SJed Brown PetscFunctionReturn(0); 7729a3665c6SJed Brown } 7739a3665c6SJed Brown 7740bad9183SKris Buschelman /*MC 775fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 7760bad9183SKris Buschelman intended for use constructing orderings and partitionings. 7770bad9183SKris Buschelman 7780bad9183SKris Buschelman Level: beginner 7790bad9183SKris Buschelman 780db781477SPatrick Sanan .seealso: `MatCreateMPIAdj` 7810bad9183SKris Buschelman M*/ 7820bad9183SKris Buschelman 7838cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 784273d9f13SBarry Smith { 785273d9f13SBarry Smith Mat_MPIAdj *b; 786273d9f13SBarry Smith 787273d9f13SBarry Smith PetscFunctionBegin; 7889566063dSJacob Faibussowitsch PetscCall(PetscNewLog(B,&b)); 789b0a32e0cSBarry Smith B->data = (void*)b; 7909566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps))); 791273d9f13SBarry Smith B->assembled = PETSC_FALSE; 792273d9f13SBarry Smith 7939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj)); 7949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj)); 7959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj)); 7969566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ)); 797273d9f13SBarry Smith PetscFunctionReturn(0); 798273d9f13SBarry Smith } 799273d9f13SBarry Smith 800a23d5eceSKris Buschelman /*@C 801b44e7856SBarry Smith MatMPIAdjToSeq - Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential preconditioners) 802b44e7856SBarry Smith 803d083f849SBarry Smith Logically Collective 804b44e7856SBarry Smith 805b44e7856SBarry Smith Input Parameter: 806b44e7856SBarry Smith . A - the matrix 807b44e7856SBarry Smith 808b44e7856SBarry Smith Output Parameter: 809b44e7856SBarry Smith . B - the same matrix on all processes 810b44e7856SBarry Smith 811b44e7856SBarry Smith Level: intermediate 812b44e7856SBarry Smith 813db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()` 814b44e7856SBarry Smith @*/ 815b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq(Mat A,Mat *B) 816b44e7856SBarry Smith { 817b44e7856SBarry Smith PetscFunctionBegin; 818cac4c232SBarry Smith PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B)); 819b44e7856SBarry Smith PetscFunctionReturn(0); 820b44e7856SBarry Smith } 821b44e7856SBarry Smith 822b44e7856SBarry Smith /*@C 823a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 824a23d5eceSKris Buschelman 825d083f849SBarry Smith Logically Collective 826a23d5eceSKris Buschelman 827a23d5eceSKris Buschelman Input Parameters: 828a23d5eceSKris Buschelman + A - the matrix 829a23d5eceSKris Buschelman . i - the indices into j for the start of each row 830a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 831a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 832a23d5eceSKris Buschelman - values - [optional] edge weights 833a23d5eceSKris Buschelman 834a23d5eceSKris Buschelman Level: intermediate 835a23d5eceSKris Buschelman 836db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()` 837a23d5eceSKris Buschelman @*/ 8387087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 839273d9f13SBarry Smith { 840273d9f13SBarry Smith PetscFunctionBegin; 841cac4c232SBarry Smith PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values)); 842273d9f13SBarry Smith PetscFunctionReturn(0); 843273d9f13SBarry Smith } 844273d9f13SBarry Smith 845b97cf49bSBarry Smith /*@C 8463eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 847b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 848b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 849b97cf49bSBarry Smith 850d083f849SBarry Smith Collective 851ef5ee4d1SLois Curfman McInnes 852b97cf49bSBarry Smith Input Parameters: 853c2e958feSBarry Smith + comm - MPI communicator 8540752156aSBarry Smith . m - number of local rows 85558ec18a5SLisandro Dalcin . N - number of global columns 856b97cf49bSBarry Smith . i - the indices into j for the start of each row 857329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 858ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 859329f5518SBarry Smith - values -[optional] edge weights 860b97cf49bSBarry Smith 861b97cf49bSBarry Smith Output Parameter: 862b97cf49bSBarry Smith . A - the matrix 863b97cf49bSBarry Smith 86415091d37SBarry Smith Level: intermediate 86515091d37SBarry Smith 86695452b02SPatrick Sanan Notes: 86795452b02SPatrick Sanan This matrix object does not support most matrix operations, include 8684bc6d8c0SBarry Smith MatSetValues(). 869c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 870ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 8711198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 872327e1a73SBarry Smith Should not include the matrix diagonals. 873b97cf49bSBarry Smith 874e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 875ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 876ededeb1aSvictorle 877ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 878b97cf49bSBarry Smith 879db781477SPatrick Sanan .seealso: `MatCreate()`, `MatConvert()`, `MatGetOrdering()` 880b97cf49bSBarry Smith @*/ 8817087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 882b97cf49bSBarry Smith { 883433994e6SBarry Smith PetscFunctionBegin; 8849566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,A)); 8859566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N)); 8869566063dSJacob Faibussowitsch PetscCall(MatSetType(*A,MATMPIADJ)); 8879566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(*A,i,j,values)); 8883a40ed3dSBarry Smith PetscFunctionReturn(0); 889b97cf49bSBarry Smith } 890