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 */ 155252a1336SBarry Smith // if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscMalloc1(n,submat)); 1566a7b62d2SBarry Smith 15740244768SBarry Smith for (i=0; i<n; i++) { 15840244768SBarry Smith if (subcomm) { 1599566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)irow[i],&scomm_row)); 1609566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)icol[i],&scomm_col)); 1619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row,scomm_col,&issame)); 16208401ef6SPierre 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"); 1639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame)); 16408401ef6SPierre 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"); 16540244768SBarry Smith } else { 16640244768SBarry Smith scomm_row = PETSC_COMM_SELF; 16740244768SBarry Smith } 16840244768SBarry Smith /*get sub-matrix data*/ 169f4259b30SLisandro Dalcin sxadj=NULL; sadjncy=NULL; svalues=NULL; 1709566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues)); 1719566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i],&irow_n)); 1729566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i],&icol_n)); 1739566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irow[i],&irow_indices)); 1749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indices,irow_indices,irow_n)); 1759566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i],&irow_indices)); 1769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icol[i],&icol_indices)); 1779566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indices+irow_n,icol_indices,icol_n)); 1789566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i],&icol_indices)); 17940244768SBarry Smith nindx = irow_n+icol_n; 1809566063dSJacob Faibussowitsch PetscCall(PetscSortRemoveDupsInt(&nindx,indices)); 18140244768SBarry Smith /* renumber columns */ 18240244768SBarry Smith for (j=0; j<irow_n; j++) { 18340244768SBarry Smith for (k=sxadj[j]; k<sxadj[j+1]; k++) { 1849566063dSJacob Faibussowitsch PetscCall(PetscFindInt(sadjncy[k],nindx,indices,&loc)); 18508401ef6SPierre Jolivet PetscCheck(loc>=0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %" PetscInt_FMT,sadjncy[k]); 18640244768SBarry Smith sadjncy[k] = loc; 18740244768SBarry Smith } 18840244768SBarry Smith } 18940244768SBarry Smith if (scall==MAT_INITIAL_MATRIX) { 1909566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i])); 19140244768SBarry Smith } else { 19240244768SBarry Smith Mat sadj = *(submat[i]); 19340244768SBarry Smith Mat_MPIAdj *sa = (Mat_MPIAdj*)((sadj)->data); 1949566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sadj,&scomm_mat)); 1959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row,scomm_mat,&issame)); 19608401ef6SPierre Jolivet PetscCheck(issame == MPI_IDENT,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"submatrix must have the same comm as the col index set"); 1979566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sa->i,sxadj,irow_n+1)); 1989566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sa->j,sadjncy,sxadj[irow_n])); 1999566063dSJacob Faibussowitsch if (svalues) PetscCall(PetscArraycpy(sa->values,svalues,sxadj[irow_n])); 2009566063dSJacob Faibussowitsch PetscCall(PetscFree(sxadj)); 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(sadjncy)); 202252a1336SBarry Smith PetscCall(PetscFree(svalues)); 20340244768SBarry Smith } 20440244768SBarry Smith } 2059566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 20640244768SBarry Smith PetscFunctionReturn(0); 20740244768SBarry Smith } 20840244768SBarry Smith 2097dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatricesMPI_MPIAdj(Mat mat,PetscInt n, const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 21040244768SBarry Smith { 21140244768SBarry Smith /*get sub-matrices across a sub communicator */ 21240244768SBarry Smith PetscFunctionBegin; 2139566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat)); 21440244768SBarry Smith PetscFunctionReturn(0); 21540244768SBarry Smith } 21640244768SBarry Smith 2177dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 21840244768SBarry Smith { 21940244768SBarry Smith PetscFunctionBegin; 22040244768SBarry Smith /*get sub-matrices based on PETSC_COMM_SELF */ 2219566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat)); 22240244768SBarry Smith PetscFunctionReturn(0); 22340244768SBarry Smith } 22440244768SBarry Smith 22540244768SBarry Smith static PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 226b97cf49bSBarry Smith { 2273eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 228d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 2292dcb1b2aSMatthew Knepley const char *name; 230ce1411ecSBarry Smith PetscViewerFormat format; 231b97cf49bSBarry Smith 232433994e6SBarry Smith PetscFunctionBegin; 2339566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)A,&name)); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 235fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 2363a40ed3dSBarry Smith PetscFunctionReturn(0); 237f7d195e4SLawrence Mitchell } else { 238f7d195e4SLawrence Mitchell PetscCheck(format != PETSC_VIEWER_ASCII_MATLAB,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported"); 2399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 241b97cf49bSBarry Smith for (i=0; i<m; i++) { 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"row %" PetscInt_FMT ":",i+A->rmap->rstart)); 243b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 2440170919eSFande Kong if (a->values) { 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," (%" PetscInt_FMT ", %" PetscInt_FMT ") ",a->j[j], a->values[j])); 2460170919eSFande Kong } else { 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT " ",a->j[j])); 2480752156aSBarry Smith } 2490170919eSFande Kong } 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n")); 251b97cf49bSBarry Smith } 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2557b23a99aSBarry Smith } 2563a40ed3dSBarry Smith PetscFunctionReturn(0); 257b97cf49bSBarry Smith } 258b97cf49bSBarry Smith 25940244768SBarry Smith static PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 260b97cf49bSBarry Smith { 261ace3abfcSBarry Smith PetscBool iascii; 262b97cf49bSBarry Smith 263433994e6SBarry Smith PetscFunctionBegin; 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 2651baa6e33SBarry Smith if (iascii) PetscCall(MatView_MPIAdj_ASCII(A,viewer)); 2663a40ed3dSBarry Smith PetscFunctionReturn(0); 267b97cf49bSBarry Smith } 268b97cf49bSBarry Smith 26940244768SBarry Smith static PetscErrorCode MatDestroy_MPIAdj(Mat mat) 270b97cf49bSBarry Smith { 2713eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 27294d884c6SBarry Smith 27394d884c6SBarry Smith PetscFunctionBegin; 274aa482453SBarry Smith #if defined(PETSC_USE_LOG) 275c0aa6a63SJacob Faibussowitsch PetscLogObjectState((PetscObject)mat,"Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT,mat->rmap->n,mat->cmap->n,a->nz); 276b97cf49bSBarry Smith #endif 2779566063dSJacob Faibussowitsch PetscCall(PetscFree(a->diag)); 2780400557aSBarry Smith if (a->freeaij) { 27914196391SBarry Smith if (a->freeaijwithfree) { 28014196391SBarry Smith if (a->i) free(a->i); 28114196391SBarry Smith if (a->j) free(a->j); 28214196391SBarry Smith } else { 2839566063dSJacob Faibussowitsch PetscCall(PetscFree(a->i)); 2849566063dSJacob Faibussowitsch PetscCall(PetscFree(a->j)); 2859566063dSJacob Faibussowitsch PetscCall(PetscFree(a->values)); 28614196391SBarry Smith } 2870400557aSBarry Smith } 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 2899566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 2909566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)mat,NULL)); 2919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL)); 2929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL)); 2932e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjToSeq_C",NULL)); 294252a1336SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjToSeqRankZero_C",NULL)); 2953a40ed3dSBarry Smith PetscFunctionReturn(0); 296b97cf49bSBarry Smith } 297b97cf49bSBarry Smith 29840244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 299b97cf49bSBarry Smith { 3003eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 301b97cf49bSBarry Smith 302433994e6SBarry Smith PetscFunctionBegin; 30312c028f9SKris Buschelman switch (op) { 30477e54ba9SKris Buschelman case MAT_SYMMETRIC: 30512c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3069a4540c5SBarry Smith case MAT_HERMITIAN: 307b94d7dedSBarry Smith case MAT_SPD: 3084e0d8c25SBarry Smith a->symmetric = flg; 3099a4540c5SBarry Smith break; 3109a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 311b94d7dedSBarry Smith case MAT_STRUCTURAL_SYMMETRY_ETERNAL: 312b94d7dedSBarry Smith case MAT_SPD_ETERNAL: 3139a4540c5SBarry Smith break; 31412c028f9SKris Buschelman default: 3159566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"Option %s ignored\n",MatOptions[op])); 31612c028f9SKris Buschelman break; 317b97cf49bSBarry Smith } 3183a40ed3dSBarry Smith PetscFunctionReturn(0); 319b97cf49bSBarry Smith } 320b97cf49bSBarry Smith 32140244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 322b97cf49bSBarry Smith { 3233eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 324b97cf49bSBarry Smith 325433994e6SBarry Smith PetscFunctionBegin; 326d0f46423SBarry Smith row -= A->rmap->rstart; 327aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 328b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 3292b1d8763SJed Brown if (v) { 3302b1d8763SJed Brown PetscInt j; 3312b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 3329566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 3332b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz); 3349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->rowvalues_alloc,&a->rowvalues)); 335b97cf49bSBarry Smith } 33691f8cafeSFande Kong for (j=0; j<*nz; j++) { 33791f8cafeSFande Kong a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0; 33891f8cafeSFande Kong } 3392b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 340b97cf49bSBarry Smith } 34192b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 3423a40ed3dSBarry Smith PetscFunctionReturn(0); 343b97cf49bSBarry Smith } 344b97cf49bSBarry Smith 34540244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 346b97cf49bSBarry Smith { 347433994e6SBarry Smith PetscFunctionBegin; 3483a40ed3dSBarry Smith PetscFunctionReturn(0); 349b97cf49bSBarry Smith } 350b97cf49bSBarry Smith 35140244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 352b97cf49bSBarry Smith { 3533eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 354ace3abfcSBarry Smith PetscBool flag; 355b97cf49bSBarry Smith 356433994e6SBarry Smith PetscFunctionBegin; 357b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 358d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 3590f5bd95cSBarry Smith flag = PETSC_FALSE; 360b97cf49bSBarry Smith } 361b97cf49bSBarry Smith 362b97cf49bSBarry Smith /* if the a->i are the same */ 3639566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(a->i,b->i,A->rmap->n+1,&flag)); 364b97cf49bSBarry Smith 365b97cf49bSBarry Smith /* if a->j are the same */ 3669566063dSJacob Faibussowitsch PetscCall(PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag)); 367b97cf49bSBarry Smith 3681c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 3693a40ed3dSBarry Smith PetscFunctionReturn(0); 370b97cf49bSBarry Smith } 371b97cf49bSBarry Smith 37240244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 3736cd91f64SBarry Smith { 374b24ad042SBarry Smith PetscInt i; 3756cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3761a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 3776cd91f64SBarry Smith 3786cd91f64SBarry Smith PetscFunctionBegin; 379d0f46423SBarry Smith *m = A->rmap->n; 3806cd91f64SBarry Smith *ia = a->i; 3816cd91f64SBarry Smith *ja = a->j; 3826cd91f64SBarry Smith *done = PETSC_TRUE; 383d42735a1SBarry Smith if (oshift) { 384d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 385d42735a1SBarry Smith (*ja)[i]++; 386d42735a1SBarry Smith } 387d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 388d42735a1SBarry Smith } 389d42735a1SBarry Smith PetscFunctionReturn(0); 390d42735a1SBarry Smith } 391d42735a1SBarry Smith 39240244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 393d42735a1SBarry Smith { 394b24ad042SBarry Smith PetscInt i; 395d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3961a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 397d42735a1SBarry Smith 398d42735a1SBarry Smith PetscFunctionBegin; 399aed4548fSBarry Smith PetscCheck(!ia || a->i == *ia,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 400aed4548fSBarry Smith PetscCheck(!ja || a->j == *ja,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 401d42735a1SBarry Smith if (oshift) { 40228b400f6SJacob Faibussowitsch PetscCheck(ia,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument"); 40328b400f6SJacob Faibussowitsch PetscCheck(ja,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument"); 404d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 405d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 406d42735a1SBarry Smith (*ja)[i]--; 407d42735a1SBarry Smith } 408d42735a1SBarry Smith } 4096cd91f64SBarry Smith PetscFunctionReturn(0); 4106cd91f64SBarry Smith } 411b97cf49bSBarry Smith 41219fd82e9SBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 41317667f90SBarry Smith { 41417667f90SBarry Smith Mat B; 41517667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 41617667f90SBarry Smith const PetscInt *rj; 41717667f90SBarry Smith const PetscScalar *ra; 41817667f90SBarry Smith MPI_Comm comm; 41917667f90SBarry Smith 42017667f90SBarry Smith PetscFunctionBegin; 4219566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 4229566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 4239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,NULL)); 42417667f90SBarry Smith 42517667f90SBarry Smith /* count the number of nonzeros per row */ 42617667f90SBarry Smith for (i=0; i<m; i++) { 4279566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i+rstart,&len,&rj,NULL)); 4285ee9ba1cSJed Brown for (j=0; j<len; j++) { 4295ee9ba1cSJed Brown if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 4305ee9ba1cSJed Brown } 43117667f90SBarry Smith nzeros += len; 4329566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i+rstart,&len,&rj,NULL)); 43317667f90SBarry Smith } 43417667f90SBarry Smith 43517667f90SBarry Smith /* malloc space for nonzeros */ 4369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros+1,&a)); 4379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N+1,&ia)); 4389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros+1,&ja)); 43917667f90SBarry Smith 44017667f90SBarry Smith nzeros = 0; 44117667f90SBarry Smith ia[0] = 0; 44217667f90SBarry Smith for (i=0; i<m; i++) { 4439566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i+rstart,&len,&rj,&ra)); 44417667f90SBarry Smith cnt = 0; 44517667f90SBarry Smith for (j=0; j<len; j++) { 4465ee9ba1cSJed Brown if (rj[j] != i+rstart) { /* if not diagonal */ 44717667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 44817667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 44917667f90SBarry Smith } 4505ee9ba1cSJed Brown } 4519566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i+rstart,&len,&rj,&ra)); 45217667f90SBarry Smith nzeros += cnt; 45317667f90SBarry Smith ia[i+1] = nzeros; 45417667f90SBarry Smith } 45517667f90SBarry Smith 4569566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A,&comm)); 4579566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&B)); 4589566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N)); 4599566063dSJacob Faibussowitsch PetscCall(MatSetType(B,type)); 4609566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(B,ia,ja,a)); 46117667f90SBarry Smith 462511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 4639566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A,&B)); 46417667f90SBarry Smith } else { 46517667f90SBarry Smith *newmat = B; 46617667f90SBarry Smith } 46717667f90SBarry Smith PetscFunctionReturn(0); 46817667f90SBarry Smith } 46917667f90SBarry Smith 4706a09307cSBarry Smith PetscErrorCode MatSetValues_MPIAdj(Mat A,PetscInt m,const PetscInt *rows,PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode im) 4716a09307cSBarry Smith { 4726a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 4736a09307cSBarry Smith PetscInt rStart, rEnd, cStart, cEnd; 4746a09307cSBarry Smith 4756a09307cSBarry Smith PetscFunctionBegin; 4766a09307cSBarry Smith PetscCheck(!adj->i,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is already assembled, cannot change its values"); 4776a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 4786a09307cSBarry Smith PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd)); 4796a09307cSBarry Smith if (!adj->ht) { 4806a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht)); 4816a09307cSBarry Smith PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject) A), 1, &A->stash)); 4826a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 4836a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 4846a09307cSBarry Smith } 4856a09307cSBarry Smith for (PetscInt r = 0; r < m; ++r) { 4866a09307cSBarry Smith PetscHashIJKey key; 4876a09307cSBarry Smith 4886a09307cSBarry Smith key.i = rows[r]; 4896a09307cSBarry Smith if (key.i < 0) continue; 4906a09307cSBarry Smith if ((key.i < rStart) || (key.i >= rEnd)) { 4916a09307cSBarry Smith PetscCall(MatStashValuesRow_Private(&A->stash, key.i, n, cols, values, PETSC_FALSE)); 4926a09307cSBarry Smith } else { 4936a09307cSBarry Smith for (PetscInt c = 0; c < n; ++c) { 4946a09307cSBarry Smith key.j = cols[c]; 4956a09307cSBarry Smith if (key.j < 0 || key.i == key.j) continue; 4966a09307cSBarry Smith PetscCall(PetscHSetIJAdd(adj->ht, key)); 4976a09307cSBarry Smith } 4986a09307cSBarry Smith } 4996a09307cSBarry Smith } 5006a09307cSBarry Smith PetscFunctionReturn(0); 5016a09307cSBarry Smith } 5026a09307cSBarry Smith 5036a09307cSBarry Smith PetscErrorCode MatAssemblyBegin_MPIAdj(Mat A, MatAssemblyType type) 5046a09307cSBarry Smith { 5056a09307cSBarry Smith PetscInt nstash, reallocs; 5066a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 5076a09307cSBarry Smith 5086a09307cSBarry Smith PetscFunctionBegin; 5096a09307cSBarry Smith if (!adj->ht) { 5106a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht)); 5116a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 5126a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 5136a09307cSBarry Smith } 5146a09307cSBarry Smith PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range)); 5156a09307cSBarry Smith PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs)); 5166a09307cSBarry Smith PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs)); 5176a09307cSBarry Smith PetscFunctionReturn(0); 5186a09307cSBarry Smith } 5196a09307cSBarry Smith 5206a09307cSBarry Smith PetscErrorCode MatAssemblyEnd_MPIAdj(Mat A, MatAssemblyType type) 5216a09307cSBarry Smith { 5226a09307cSBarry Smith PetscScalar *val; 5236a09307cSBarry Smith PetscInt *row, *col,m,rstart,*rowstarts; 5246a09307cSBarry Smith PetscInt i, j, ncols, flg, nz; 5256a09307cSBarry Smith PetscMPIInt n; 5266a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 5276a09307cSBarry Smith PetscHashIter hi; 5286a09307cSBarry Smith PetscHashIJKey key; 5296a09307cSBarry Smith PetscHSetIJ ht = adj->ht; 5306a09307cSBarry Smith 5316a09307cSBarry Smith PetscFunctionBegin; 5326a09307cSBarry Smith while (1) { 5336a09307cSBarry Smith PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg)); 5346a09307cSBarry Smith if (!flg) break; 5356a09307cSBarry Smith 5366a09307cSBarry Smith for (i = 0; i < n;) { 5376a09307cSBarry Smith /* Identify the consecutive vals belonging to the same row */ 5386a09307cSBarry Smith for (j = i, rstart = row[j]; j < n; j++) { 5396a09307cSBarry Smith if (row[j] != rstart) break; 5406a09307cSBarry Smith } 5416a09307cSBarry Smith if (j < n) ncols = j-i; 5426a09307cSBarry Smith else ncols = n-i; 5436a09307cSBarry Smith /* Set all these values with a single function call */ 5446a09307cSBarry Smith PetscCall(MatSetValues_MPIAdj(A, 1, row+i, ncols, col+i, val+i, INSERT_VALUES)); 5456a09307cSBarry Smith i = j; 5466a09307cSBarry Smith } 5476a09307cSBarry Smith } 5486a09307cSBarry Smith PetscCall(MatStashScatterEnd_Private(&A->stash)); 5496a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&A->stash)); 5506a09307cSBarry Smith 5516a09307cSBarry Smith PetscCall(MatGetLocalSize(A,&m,NULL)); 5526a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A,&rstart,NULL)); 5536a09307cSBarry Smith PetscCall(PetscCalloc1(m+1,&rowstarts)); 5546a09307cSBarry Smith PetscHashIterBegin(ht,hi); 5556a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht,hi);) { 5566a09307cSBarry Smith PetscHashIterGetKey(ht,hi,key); 5576a09307cSBarry Smith rowstarts[key.i - rstart + 1]++; 5586a09307cSBarry Smith PetscHashIterNext(ht,hi); 5596a09307cSBarry Smith } 5606a09307cSBarry Smith for (i=1; i<m+1; i++) { 5616a09307cSBarry Smith rowstarts[i] = rowstarts[i-1] + rowstarts[i]; 5626a09307cSBarry Smith } 5636a09307cSBarry Smith 5646a09307cSBarry Smith PetscCall(PetscHSetIJGetSize(ht,&nz)); 5656a09307cSBarry Smith PetscCall(PetscMalloc1(nz,&col)); 5666a09307cSBarry Smith PetscHashIterBegin(ht,hi); 5676a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht,hi);) { 5686a09307cSBarry Smith PetscHashIterGetKey(ht,hi,key); 5696a09307cSBarry Smith col[rowstarts[key.i - rstart]++] = key.j; 5706a09307cSBarry Smith PetscHashIterNext(ht,hi); 5716a09307cSBarry Smith } 5726a09307cSBarry Smith PetscCall(PetscHSetIJDestroy(&ht)); 5736a09307cSBarry Smith for (i=m; i>0; i--) { 5746a09307cSBarry Smith rowstarts[i] = rowstarts[i-1]; 5756a09307cSBarry Smith } 5766a09307cSBarry Smith rowstarts[0] = 0; 5776a09307cSBarry Smith 5786a09307cSBarry Smith for (PetscInt i=0; i<m; i++) { 5796a09307cSBarry Smith PetscCall(PetscSortInt(rowstarts[i+1]-rowstarts[i],&col[rowstarts[i]])); 5806a09307cSBarry Smith } 5816a09307cSBarry Smith 5826a09307cSBarry Smith adj->i = rowstarts; 5836a09307cSBarry Smith adj->j = col; 584252a1336SBarry Smith adj->nz = rowstarts[m]; 5856a09307cSBarry Smith adj->freeaij = PETSC_TRUE; 5866a09307cSBarry Smith PetscFunctionReturn(0); 5876a09307cSBarry Smith } 5886a09307cSBarry Smith 589b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 5906a09307cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIAdj, 5913eda8832SBarry Smith MatGetRow_MPIAdj, 5923eda8832SBarry Smith MatRestoreRow_MPIAdj, 593f4259b30SLisandro Dalcin NULL, 594f4259b30SLisandro Dalcin /* 4*/ NULL, 595f4259b30SLisandro Dalcin NULL, 596f4259b30SLisandro Dalcin NULL, 597f4259b30SLisandro Dalcin NULL, 598f4259b30SLisandro Dalcin NULL, 599f4259b30SLisandro Dalcin NULL, 600f4259b30SLisandro Dalcin /*10*/ NULL, 601f4259b30SLisandro Dalcin NULL, 602f4259b30SLisandro Dalcin NULL, 603f4259b30SLisandro Dalcin NULL, 604f4259b30SLisandro Dalcin NULL, 605f4259b30SLisandro Dalcin /*15*/ NULL, 6063eda8832SBarry Smith MatEqual_MPIAdj, 607f4259b30SLisandro Dalcin NULL, 608f4259b30SLisandro Dalcin NULL, 609f4259b30SLisandro Dalcin NULL, 6106a09307cSBarry Smith /*20*/ MatAssemblyBegin_MPIAdj, 6116a09307cSBarry Smith MatAssemblyEnd_MPIAdj, 6123eda8832SBarry Smith MatSetOption_MPIAdj, 613f4259b30SLisandro Dalcin NULL, 614f4259b30SLisandro Dalcin /*24*/ NULL, 615f4259b30SLisandro Dalcin NULL, 616f4259b30SLisandro Dalcin NULL, 617f4259b30SLisandro Dalcin NULL, 618f4259b30SLisandro Dalcin NULL, 619f4259b30SLisandro Dalcin /*29*/ NULL, 620f4259b30SLisandro Dalcin NULL, 621f4259b30SLisandro Dalcin NULL, 622f4259b30SLisandro Dalcin NULL, 623f4259b30SLisandro Dalcin NULL, 624f4259b30SLisandro Dalcin /*34*/ NULL, 625f4259b30SLisandro Dalcin NULL, 626f4259b30SLisandro Dalcin NULL, 627f4259b30SLisandro Dalcin NULL, 628f4259b30SLisandro Dalcin NULL, 629f4259b30SLisandro Dalcin /*39*/ NULL, 6307dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 631f4259b30SLisandro Dalcin NULL, 632f4259b30SLisandro Dalcin NULL, 633f4259b30SLisandro Dalcin NULL, 634f4259b30SLisandro Dalcin /*44*/ NULL, 635f4259b30SLisandro Dalcin NULL, 6367d68702bSBarry Smith MatShift_Basic, 637f4259b30SLisandro Dalcin NULL, 638f4259b30SLisandro Dalcin NULL, 639f4259b30SLisandro Dalcin /*49*/ NULL, 6406cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 641d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 642f4259b30SLisandro Dalcin NULL, 643f4259b30SLisandro Dalcin NULL, 644f4259b30SLisandro Dalcin /*54*/ NULL, 645f4259b30SLisandro Dalcin NULL, 646f4259b30SLisandro Dalcin NULL, 647f4259b30SLisandro Dalcin NULL, 648f4259b30SLisandro Dalcin NULL, 649f4259b30SLisandro Dalcin /*59*/ NULL, 650b9b97703SBarry Smith MatDestroy_MPIAdj, 651b9b97703SBarry Smith MatView_MPIAdj, 65217667f90SBarry Smith MatConvertFrom_MPIAdj, 653f4259b30SLisandro Dalcin NULL, 654f4259b30SLisandro Dalcin /*64*/ NULL, 655f4259b30SLisandro Dalcin NULL, 656f4259b30SLisandro Dalcin NULL, 657f4259b30SLisandro Dalcin NULL, 658f4259b30SLisandro Dalcin NULL, 659f4259b30SLisandro Dalcin /*69*/ NULL, 660f4259b30SLisandro Dalcin NULL, 661f4259b30SLisandro Dalcin NULL, 662f4259b30SLisandro Dalcin NULL, 663f4259b30SLisandro Dalcin NULL, 664f4259b30SLisandro Dalcin /*74*/ NULL, 665f4259b30SLisandro Dalcin NULL, 666f4259b30SLisandro Dalcin NULL, 667f4259b30SLisandro Dalcin NULL, 668f4259b30SLisandro Dalcin NULL, 669f4259b30SLisandro Dalcin /*79*/ NULL, 670f4259b30SLisandro Dalcin NULL, 671f4259b30SLisandro Dalcin NULL, 672f4259b30SLisandro Dalcin NULL, 673f4259b30SLisandro Dalcin NULL, 674f4259b30SLisandro Dalcin /*84*/ NULL, 675f4259b30SLisandro Dalcin NULL, 676f4259b30SLisandro Dalcin NULL, 677f4259b30SLisandro Dalcin NULL, 678f4259b30SLisandro Dalcin NULL, 679f4259b30SLisandro Dalcin /*89*/ NULL, 680f4259b30SLisandro Dalcin NULL, 681f4259b30SLisandro Dalcin NULL, 682f4259b30SLisandro Dalcin NULL, 683f4259b30SLisandro Dalcin NULL, 684f4259b30SLisandro Dalcin /*94*/ NULL, 685f4259b30SLisandro Dalcin NULL, 686f4259b30SLisandro Dalcin NULL, 687f4259b30SLisandro Dalcin NULL, 688f4259b30SLisandro Dalcin NULL, 689f4259b30SLisandro Dalcin /*99*/ NULL, 690f4259b30SLisandro Dalcin NULL, 691f4259b30SLisandro Dalcin NULL, 692f4259b30SLisandro Dalcin NULL, 693f4259b30SLisandro Dalcin NULL, 694f4259b30SLisandro Dalcin /*104*/ NULL, 695f4259b30SLisandro Dalcin NULL, 696f4259b30SLisandro Dalcin NULL, 697f4259b30SLisandro Dalcin NULL, 698f4259b30SLisandro Dalcin NULL, 699f4259b30SLisandro Dalcin /*109*/ NULL, 700f4259b30SLisandro Dalcin NULL, 701f4259b30SLisandro Dalcin NULL, 702f4259b30SLisandro Dalcin NULL, 703f4259b30SLisandro Dalcin NULL, 704f4259b30SLisandro Dalcin /*114*/ NULL, 705f4259b30SLisandro Dalcin NULL, 706f4259b30SLisandro Dalcin NULL, 707f4259b30SLisandro Dalcin NULL, 708f4259b30SLisandro Dalcin NULL, 709f4259b30SLisandro Dalcin /*119*/ NULL, 710f4259b30SLisandro Dalcin NULL, 711f4259b30SLisandro Dalcin NULL, 712f4259b30SLisandro Dalcin NULL, 713f4259b30SLisandro Dalcin NULL, 714f4259b30SLisandro Dalcin /*124*/ NULL, 715f4259b30SLisandro Dalcin NULL, 716f4259b30SLisandro Dalcin NULL, 717f4259b30SLisandro Dalcin NULL, 7187dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 719f4259b30SLisandro Dalcin /*129*/ NULL, 720f4259b30SLisandro Dalcin NULL, 721f4259b30SLisandro Dalcin NULL, 722f4259b30SLisandro Dalcin NULL, 723f4259b30SLisandro Dalcin NULL, 724f4259b30SLisandro Dalcin /*134*/ NULL, 725f4259b30SLisandro Dalcin NULL, 726f4259b30SLisandro Dalcin NULL, 727f4259b30SLisandro Dalcin NULL, 728f4259b30SLisandro Dalcin NULL, 729f4259b30SLisandro Dalcin /*139*/ NULL, 730f4259b30SLisandro Dalcin NULL, 731d70f29a3SPierre Jolivet NULL, 732d70f29a3SPierre Jolivet NULL, 733d70f29a3SPierre Jolivet NULL, 734d70f29a3SPierre Jolivet /*144*/ NULL, 735d70f29a3SPierre Jolivet NULL, 736d70f29a3SPierre Jolivet NULL, 73799a7f59eSMark Adams NULL, 73899a7f59eSMark Adams NULL, 739*7fb60732SBarry Smith NULL, 740*7fb60732SBarry Smith /*150*/ NULL 7413964eb88SJed Brown }; 742b97cf49bSBarry Smith 743f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 744a23d5eceSKris Buschelman { 745a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 746e895ccc0SFande Kong PetscBool useedgeweights; 747a23d5eceSKris Buschelman 748a23d5eceSKris Buschelman PetscFunctionBegin; 7499566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->rmap)); 7509566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->cmap)); 751e895ccc0SFande Kong if (values) useedgeweights = PETSC_TRUE; else useedgeweights = PETSC_FALSE; 752e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */ 7531c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((int*)&useedgeweights,(int*)&b->useedgeweights,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)B))); 75458ec18a5SLisandro Dalcin 75576bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 75676bd3646SJed Brown PetscInt ii; 75776bd3646SJed Brown 75808401ef6SPierre 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]); 759d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 760aed4548fSBarry 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]); 761a23d5eceSKris Buschelman } 762d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 763aed4548fSBarry 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]); 764a23d5eceSKris Buschelman } 76576bd3646SJed Brown } 766a23d5eceSKris Buschelman b->j = j; 767a23d5eceSKris Buschelman b->i = i; 768a23d5eceSKris Buschelman b->values = values; 769a23d5eceSKris Buschelman 770d0f46423SBarry Smith b->nz = i[B->rmap->n]; 771f4259b30SLisandro Dalcin b->diag = NULL; 772a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 773a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 774a23d5eceSKris Buschelman 7756a09307cSBarry Smith B->ops->assemblybegin = NULL; 7766a09307cSBarry Smith B->ops->assemblyend = NULL; 7779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 7789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 7796a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&B->stash)); 780a23d5eceSKris Buschelman PetscFunctionReturn(0); 781a23d5eceSKris Buschelman } 782b97cf49bSBarry Smith 783f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 7849a3665c6SJed Brown { 7859a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 7869a3665c6SJed Brown const PetscInt *ranges; 7879a3665c6SJed Brown MPI_Comm acomm,bcomm; 7889a3665c6SJed Brown MPI_Group agroup,bgroup; 7899a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 7909a3665c6SJed Brown 7919a3665c6SJed Brown PetscFunctionBegin; 7920298fd71SBarry Smith *B = NULL; 7939566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A,&acomm)); 7949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm,&size)); 7959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm,&rank)); 7969566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A,&ranges)); 7979a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 7989a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 7999a3665c6SJed Brown } 8009a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 8019566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 8029a3665c6SJed Brown *B = A; 8039a3665c6SJed Brown PetscFunctionReturn(0); 8049a3665c6SJed Brown } 8059a3665c6SJed Brown 8069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nranks,&ranks)); 8079a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 8089a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 8099a3665c6SJed Brown } 8109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(acomm,&agroup)); 8119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(agroup,nranks,ranks,&bgroup)); 8129566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks)); 8139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(acomm,bgroup,&bcomm)); 8149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&agroup)); 8159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&bgroup)); 8169a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 8179a3665c6SJed Brown PetscInt m,N; 8189a3665c6SJed Brown Mat_MPIAdj *b; 8199566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 8209566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 8219566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B)); 8229a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 8239a3665c6SJed Brown b->freeaij = PETSC_FALSE; 8249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&bcomm)); 8259a3665c6SJed Brown } 8269a3665c6SJed Brown PetscFunctionReturn(0); 8279a3665c6SJed Brown } 8289a3665c6SJed Brown 829b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B) 830b44e7856SBarry Smith { 831b44e7856SBarry Smith PetscInt M,N,*II,*J,NZ,nz,m,nzstart,i; 832b44e7856SBarry Smith PetscInt *Values = NULL; 833b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 834b44e7856SBarry Smith PetscMPIInt mnz,mm,*allnz,*allm,size,*dispnz,*dispm; 835b44e7856SBarry Smith 836b44e7856SBarry Smith PetscFunctionBegin; 8379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 8389566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 8399566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 840b44e7856SBarry Smith nz = adj->nz; 84108401ef6SPierre 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]); 8429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A))); 843d4e69552SBarry Smith 8449566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nz,&mnz)); 8459566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size,&allnz,size,&dispnz)); 8469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A))); 847b44e7856SBarry Smith dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1]; 848b44e7856SBarry Smith if (adj->values) { 8499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ,&Values)); 8509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A))); 851b44e7856SBarry Smith } 8529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ,&J)); 8539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A))); 8549566063dSJacob Faibussowitsch PetscCall(PetscFree2(allnz,dispnz)); 8559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A))); 856b44e7856SBarry Smith nzstart -= nz; 857b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 858b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] += nzstart; 8599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M+1,&II)); 8609566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(m,&mm)); 8619566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size,&allm,size,&dispm)); 8629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A))); 863b44e7856SBarry Smith dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1]; 8649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A))); 8659566063dSJacob Faibussowitsch PetscCall(PetscFree2(allm,dispm)); 866b44e7856SBarry Smith II[M] = NZ; 867b44e7856SBarry Smith /* shift the i[] values back */ 868b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] -= nzstart; 8699566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B)); 870b44e7856SBarry Smith PetscFunctionReturn(0); 871b44e7856SBarry Smith } 872b44e7856SBarry Smith 873252a1336SBarry Smith PetscErrorCode MatMPIAdjToSeqRankZero_MPIAdj(Mat A,Mat *B) 874252a1336SBarry Smith { 875252a1336SBarry Smith PetscInt M,N,*II,*J,NZ,nz,m,nzstart,i; 876252a1336SBarry Smith PetscInt *Values = NULL; 877252a1336SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 878252a1336SBarry Smith PetscMPIInt mnz,mm,*allnz = NULL,*allm,size,*dispnz,*dispm,rank; 879252a1336SBarry Smith 880252a1336SBarry Smith PetscFunctionBegin; 881252a1336SBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 882252a1336SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank)); 883252a1336SBarry Smith PetscCall(MatGetSize(A,&M,&N)); 884252a1336SBarry Smith PetscCall(MatGetLocalSize(A,&m,NULL)); 885252a1336SBarry Smith nz = adj->nz; 886252a1336SBarry Smith PetscCheck(adj->i[m] == nz,PETSC_COMM_SELF,PETSC_ERR_PLIB,"nz %" PetscInt_FMT " not correct i[m] %" PetscInt_FMT,nz,adj->i[m]); 887252a1336SBarry Smith PetscCallMPI(MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A))); 888252a1336SBarry Smith 889252a1336SBarry Smith PetscCall(PetscMPIIntCast(nz,&mnz)); 890252a1336SBarry Smith if (!rank) PetscCall(PetscMalloc2(size,&allnz,size,&dispnz)); 891252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mnz,1,MPI_INT,allnz,1,MPI_INT,0,PetscObjectComm((PetscObject)A))); 892252a1336SBarry Smith if (!rank) { 893252a1336SBarry Smith dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1]; 894252a1336SBarry Smith if (adj->values) { 895252a1336SBarry Smith PetscCall(PetscMalloc1(NZ,&Values)); 896252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,0,PetscObjectComm((PetscObject)A))); 897252a1336SBarry Smith } 898252a1336SBarry Smith PetscCall(PetscMalloc1(NZ,&J)); 899252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,0,PetscObjectComm((PetscObject)A))); 900252a1336SBarry Smith PetscCall(PetscFree2(allnz,dispnz)); 901252a1336SBarry Smith } else { 902252a1336SBarry Smith if (adj->values) PetscCallMPI(MPI_Gatherv(adj->values,mnz,MPIU_INT,NULL,NULL,NULL,MPIU_INT,0,PetscObjectComm((PetscObject)A))); 903252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->j,mnz,MPIU_INT,NULL,NULL,NULL,MPIU_INT,0,PetscObjectComm((PetscObject)A))); 904252a1336SBarry Smith } 905252a1336SBarry Smith PetscCallMPI(MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A))); 906252a1336SBarry Smith nzstart -= nz; 907252a1336SBarry Smith /* shift the i[] values so they will be correct after being received */ 908252a1336SBarry Smith for (i=0; i<m; i++) adj->i[i] += nzstart; 909252a1336SBarry Smith PetscCall(PetscMPIIntCast(m,&mm)); 910252a1336SBarry Smith if (!rank) { 911252a1336SBarry Smith PetscCall(PetscMalloc1(M+1,&II)); 912252a1336SBarry Smith PetscCall(PetscMalloc2(size,&allm,size,&dispm)); 913252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mm,1,MPI_INT,allm,1,MPI_INT,0,PetscObjectComm((PetscObject)A))); 914252a1336SBarry Smith dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1]; 915252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,0,PetscObjectComm((PetscObject)A))); 916252a1336SBarry Smith PetscCall(PetscFree2(allm,dispm)); 917252a1336SBarry Smith II[M] = NZ; 918252a1336SBarry Smith } else { 919252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mm,1,MPI_INT,NULL,1,MPI_INT,0,PetscObjectComm((PetscObject)A))); 920252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->i,mm,MPIU_INT,NULL,NULL,NULL,MPIU_INT,0,PetscObjectComm((PetscObject)A))); 921252a1336SBarry Smith } 922252a1336SBarry Smith /* shift the i[] values back */ 923252a1336SBarry Smith for (i=0; i<m; i++) adj->i[i] -= nzstart; 924252a1336SBarry Smith if (!rank) PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B)); 925252a1336SBarry Smith PetscFunctionReturn(0); 926252a1336SBarry Smith } 927252a1336SBarry Smith 9289a3665c6SJed Brown /*@ 9299a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 9309a3665c6SJed Brown 9319a3665c6SJed Brown Collective 9329a3665c6SJed Brown 9334165533cSJose E. Roman Input Parameter: 9349a3665c6SJed Brown . A - original MPIAdj matrix 9359a3665c6SJed Brown 9364165533cSJose E. Roman Output Parameter: 9370298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 9389a3665c6SJed Brown 9399a3665c6SJed Brown Level: developer 9409a3665c6SJed Brown 9419a3665c6SJed Brown Note: 9429a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 9439a3665c6SJed Brown 9449a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 9459a3665c6SJed Brown 946db781477SPatrick Sanan .seealso: `MatCreateMPIAdj()` 9479a3665c6SJed Brown @*/ 9489a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 9499a3665c6SJed Brown { 9509a3665c6SJed Brown PetscFunctionBegin; 9519a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 952cac4c232SBarry Smith PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B)); 9539a3665c6SJed Brown PetscFunctionReturn(0); 9549a3665c6SJed Brown } 9559a3665c6SJed Brown 9560bad9183SKris Buschelman /*MC 957fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 9580bad9183SKris Buschelman intended for use constructing orderings and partitionings. 9590bad9183SKris Buschelman 9600bad9183SKris Buschelman Level: beginner 9610bad9183SKris Buschelman 9626a09307cSBarry Smith Notes: 9636a09307cSBarry Smith You can provide values to the matrix using MatMPIAdjSetPreallocation(), MatCreateMPIAdj(), or 9646a09307cSBarry Smith by calling MatSetValues() and MatAssemblyBegin() followed by MatAssemblyEnd() 9656a09307cSBarry Smith 9666a09307cSBarry Smith .seealso: `MatCreateMPIAdj()`, `MatMPIAdjSetPreallocation()`, `MatSetValues()` 9670bad9183SKris Buschelman M*/ 9680bad9183SKris Buschelman 9698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 970273d9f13SBarry Smith { 971273d9f13SBarry Smith Mat_MPIAdj *b; 972273d9f13SBarry Smith 973273d9f13SBarry Smith PetscFunctionBegin; 9749566063dSJacob Faibussowitsch PetscCall(PetscNewLog(B,&b)); 975b0a32e0cSBarry Smith B->data = (void*)b; 9769566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps))); 977273d9f13SBarry Smith B->assembled = PETSC_FALSE; 9786a09307cSBarry Smith B->preallocated = PETSC_TRUE; /* so that MatSetValues() may be used */ 979273d9f13SBarry Smith 9809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj)); 9819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj)); 9829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj)); 983252a1336SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeqRankZero_C",MatMPIAdjToSeqRankZero_MPIAdj)); 9849566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ)); 985273d9f13SBarry Smith PetscFunctionReturn(0); 986273d9f13SBarry Smith } 987273d9f13SBarry Smith 988a23d5eceSKris Buschelman /*@C 989252a1336SBarry Smith MatMPIAdjToSeq - Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential partitioner) 990b44e7856SBarry Smith 991d083f849SBarry Smith Logically Collective 992b44e7856SBarry Smith 993b44e7856SBarry Smith Input Parameter: 994b44e7856SBarry Smith . A - the matrix 995b44e7856SBarry Smith 996b44e7856SBarry Smith Output Parameter: 997b44e7856SBarry Smith . B - the same matrix on all processes 998b44e7856SBarry Smith 999b44e7856SBarry Smith Level: intermediate 1000b44e7856SBarry Smith 1001252a1336SBarry Smith .seealso: `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeqRankZero()` 1002b44e7856SBarry Smith @*/ 1003b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq(Mat A,Mat *B) 1004b44e7856SBarry Smith { 1005b44e7856SBarry Smith PetscFunctionBegin; 1006cac4c232SBarry Smith PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B)); 1007b44e7856SBarry Smith PetscFunctionReturn(0); 1008b44e7856SBarry Smith } 1009b44e7856SBarry Smith 1010b44e7856SBarry Smith /*@C 1011252a1336SBarry Smith MatMPIAdjToSeqRankZero - Converts an parallel MPIAdj matrix to complete MPIAdj on rank zero (needed by sequential partitioner) 1012252a1336SBarry Smith 1013252a1336SBarry Smith Logically Collective 1014252a1336SBarry Smith 1015252a1336SBarry Smith Input Parameter: 1016252a1336SBarry Smith . A - the matrix 1017252a1336SBarry Smith 1018252a1336SBarry Smith Output Parameter: 1019252a1336SBarry Smith . B - the same matrix on rank zero, not set on other ranks 1020252a1336SBarry Smith 1021252a1336SBarry Smith Level: intermediate 1022252a1336SBarry Smith 1023252a1336SBarry Smith Notes: 1024252a1336SBarry Smith This routine has the advantage on systems with multiple ranks per node since only one copy of the matrix 1025252a1336SBarry Smith is stored on the first node, instead of the number of ranks copies. This can allow partitioning much larger 1026252a1336SBarry Smith paralllel graph sequentially. 1027252a1336SBarry Smith 1028252a1336SBarry Smith .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues(), MatMPIAdjToSeq() 1029252a1336SBarry Smith @*/ 1030252a1336SBarry Smith PetscErrorCode MatMPIAdjToSeqRankZero(Mat A,Mat *B) 1031252a1336SBarry Smith { 1032252a1336SBarry Smith PetscFunctionBegin; 1033252a1336SBarry Smith PetscUseMethod(A,"MatMPIAdjToSeqRankZero_C",(Mat,Mat*),(A,B)); 1034252a1336SBarry Smith PetscFunctionReturn(0); 1035252a1336SBarry Smith } 1036252a1336SBarry Smith 1037252a1336SBarry Smith /*@C 1038a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 1039a23d5eceSKris Buschelman 1040d083f849SBarry Smith Logically Collective 1041a23d5eceSKris Buschelman 1042a23d5eceSKris Buschelman Input Parameters: 1043a23d5eceSKris Buschelman + A - the matrix 1044a23d5eceSKris Buschelman . i - the indices into j for the start of each row 1045a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 1046a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 1047a23d5eceSKris Buschelman - values - [optional] edge weights 1048a23d5eceSKris Buschelman 1049a23d5eceSKris Buschelman Level: intermediate 1050a23d5eceSKris Buschelman 10516a09307cSBarry Smith .seealso: `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MATMPIADJ`, `MatCreateMPIAdj()` 1052a23d5eceSKris Buschelman @*/ 10537087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 1054273d9f13SBarry Smith { 1055273d9f13SBarry Smith PetscFunctionBegin; 1056cac4c232SBarry Smith PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values)); 1057273d9f13SBarry Smith PetscFunctionReturn(0); 1058273d9f13SBarry Smith } 1059273d9f13SBarry Smith 1060b97cf49bSBarry Smith /*@C 10613eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 1062b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 1063b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 1064b97cf49bSBarry Smith 1065d083f849SBarry Smith Collective 1066ef5ee4d1SLois Curfman McInnes 1067b97cf49bSBarry Smith Input Parameters: 1068c2e958feSBarry Smith + comm - MPI communicator 10690752156aSBarry Smith . m - number of local rows 107058ec18a5SLisandro Dalcin . N - number of global columns 1071b97cf49bSBarry Smith . i - the indices into j for the start of each row 1072329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 1073ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 1074329f5518SBarry Smith - values -[optional] edge weights 1075b97cf49bSBarry Smith 1076b97cf49bSBarry Smith Output Parameter: 1077b97cf49bSBarry Smith . A - the matrix 1078b97cf49bSBarry Smith 107915091d37SBarry Smith Level: intermediate 108015091d37SBarry Smith 108195452b02SPatrick Sanan Notes: 1082c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 1083ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 10841198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 10856a09307cSBarry Smith 10866a09307cSBarry Smith You should not include the matrix diagonals. 1087b97cf49bSBarry Smith 1088e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 10896a09307cSBarry Smith to MatConvert(), specifying a type of MATMPIADJ. 1090ededeb1aSvictorle 1091ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 1092b97cf49bSBarry Smith 10936a09307cSBarry Smith .seealso: `MatCreate()`, `MatConvert()`, `MatGetOrdering()`, `MATMPIADJ`, `MatMPIAdjSetPreallocation()` 1094b97cf49bSBarry Smith @*/ 10957087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 1096b97cf49bSBarry Smith { 1097433994e6SBarry Smith PetscFunctionBegin; 10989566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,A)); 10999566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N)); 11009566063dSJacob Faibussowitsch PetscCall(MatSetType(*A,MATMPIADJ)); 11019566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(*A,i,j,values)); 11023a40ed3dSBarry Smith PetscFunctionReturn(0); 1103b97cf49bSBarry Smith } 1104