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 */ 155*6a7b62d2SBarry Smith if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscMalloc1(n,submat)); 156*6a7b62d2SBarry 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)); 2029566063dSJacob Faibussowitsch if (svalues) 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)); 2943a40ed3dSBarry Smith PetscFunctionReturn(0); 295b97cf49bSBarry Smith } 296b97cf49bSBarry Smith 29740244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 298b97cf49bSBarry Smith { 2993eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 300b97cf49bSBarry Smith 301433994e6SBarry Smith PetscFunctionBegin; 30212c028f9SKris Buschelman switch (op) { 30377e54ba9SKris Buschelman case MAT_SYMMETRIC: 30412c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3059a4540c5SBarry Smith case MAT_HERMITIAN: 3064e0d8c25SBarry Smith a->symmetric = flg; 3079a4540c5SBarry Smith break; 3089a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 3099a4540c5SBarry Smith break; 31012c028f9SKris Buschelman default: 3119566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"Option %s ignored\n",MatOptions[op])); 31212c028f9SKris Buschelman break; 313b97cf49bSBarry Smith } 3143a40ed3dSBarry Smith PetscFunctionReturn(0); 315b97cf49bSBarry Smith } 316b97cf49bSBarry Smith 31740244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 318b97cf49bSBarry Smith { 3193eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 320b97cf49bSBarry Smith 321433994e6SBarry Smith PetscFunctionBegin; 322d0f46423SBarry Smith row -= A->rmap->rstart; 323aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 324b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 3252b1d8763SJed Brown if (v) { 3262b1d8763SJed Brown PetscInt j; 3272b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 3289566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 3292b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz); 3309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->rowvalues_alloc,&a->rowvalues)); 331b97cf49bSBarry Smith } 33291f8cafeSFande Kong for (j=0; j<*nz; j++) { 33391f8cafeSFande Kong a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0; 33491f8cafeSFande Kong } 3352b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 336b97cf49bSBarry Smith } 33792b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 3383a40ed3dSBarry Smith PetscFunctionReturn(0); 339b97cf49bSBarry Smith } 340b97cf49bSBarry Smith 34140244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 342b97cf49bSBarry Smith { 343433994e6SBarry Smith PetscFunctionBegin; 3443a40ed3dSBarry Smith PetscFunctionReturn(0); 345b97cf49bSBarry Smith } 346b97cf49bSBarry Smith 34740244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 348b97cf49bSBarry Smith { 3493eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 350ace3abfcSBarry Smith PetscBool flag; 351b97cf49bSBarry Smith 352433994e6SBarry Smith PetscFunctionBegin; 353b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 354d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 3550f5bd95cSBarry Smith flag = PETSC_FALSE; 356b97cf49bSBarry Smith } 357b97cf49bSBarry Smith 358b97cf49bSBarry Smith /* if the a->i are the same */ 3599566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(a->i,b->i,A->rmap->n+1,&flag)); 360b97cf49bSBarry Smith 361b97cf49bSBarry Smith /* if a->j are the same */ 3629566063dSJacob Faibussowitsch PetscCall(PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag)); 363b97cf49bSBarry Smith 3641c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 3653a40ed3dSBarry Smith PetscFunctionReturn(0); 366b97cf49bSBarry Smith } 367b97cf49bSBarry Smith 36840244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 3696cd91f64SBarry Smith { 370b24ad042SBarry Smith PetscInt i; 3716cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3721a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 3736cd91f64SBarry Smith 3746cd91f64SBarry Smith PetscFunctionBegin; 375d0f46423SBarry Smith *m = A->rmap->n; 3766cd91f64SBarry Smith *ia = a->i; 3776cd91f64SBarry Smith *ja = a->j; 3786cd91f64SBarry Smith *done = PETSC_TRUE; 379d42735a1SBarry Smith if (oshift) { 380d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 381d42735a1SBarry Smith (*ja)[i]++; 382d42735a1SBarry Smith } 383d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 384d42735a1SBarry Smith } 385d42735a1SBarry Smith PetscFunctionReturn(0); 386d42735a1SBarry Smith } 387d42735a1SBarry Smith 38840244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 389d42735a1SBarry Smith { 390b24ad042SBarry Smith PetscInt i; 391d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3921a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 393d42735a1SBarry Smith 394d42735a1SBarry Smith PetscFunctionBegin; 395aed4548fSBarry Smith PetscCheck(!ia || a->i == *ia,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 396aed4548fSBarry Smith PetscCheck(!ja || a->j == *ja,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 397d42735a1SBarry Smith if (oshift) { 39828b400f6SJacob Faibussowitsch PetscCheck(ia,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument"); 39928b400f6SJacob Faibussowitsch PetscCheck(ja,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument"); 400d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 401d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 402d42735a1SBarry Smith (*ja)[i]--; 403d42735a1SBarry Smith } 404d42735a1SBarry Smith } 4056cd91f64SBarry Smith PetscFunctionReturn(0); 4066cd91f64SBarry Smith } 407b97cf49bSBarry Smith 40819fd82e9SBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 40917667f90SBarry Smith { 41017667f90SBarry Smith Mat B; 41117667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 41217667f90SBarry Smith const PetscInt *rj; 41317667f90SBarry Smith const PetscScalar *ra; 41417667f90SBarry Smith MPI_Comm comm; 41517667f90SBarry Smith 41617667f90SBarry Smith PetscFunctionBegin; 4179566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 4189566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 4199566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,NULL)); 42017667f90SBarry Smith 42117667f90SBarry Smith /* count the number of nonzeros per row */ 42217667f90SBarry Smith for (i=0; i<m; i++) { 4239566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i+rstart,&len,&rj,NULL)); 4245ee9ba1cSJed Brown for (j=0; j<len; j++) { 4255ee9ba1cSJed Brown if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 4265ee9ba1cSJed Brown } 42717667f90SBarry Smith nzeros += len; 4289566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i+rstart,&len,&rj,NULL)); 42917667f90SBarry Smith } 43017667f90SBarry Smith 43117667f90SBarry Smith /* malloc space for nonzeros */ 4329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros+1,&a)); 4339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N+1,&ia)); 4349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros+1,&ja)); 43517667f90SBarry Smith 43617667f90SBarry Smith nzeros = 0; 43717667f90SBarry Smith ia[0] = 0; 43817667f90SBarry Smith for (i=0; i<m; i++) { 4399566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i+rstart,&len,&rj,&ra)); 44017667f90SBarry Smith cnt = 0; 44117667f90SBarry Smith for (j=0; j<len; j++) { 4425ee9ba1cSJed Brown if (rj[j] != i+rstart) { /* if not diagonal */ 44317667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 44417667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 44517667f90SBarry Smith } 4465ee9ba1cSJed Brown } 4479566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i+rstart,&len,&rj,&ra)); 44817667f90SBarry Smith nzeros += cnt; 44917667f90SBarry Smith ia[i+1] = nzeros; 45017667f90SBarry Smith } 45117667f90SBarry Smith 4529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A,&comm)); 4539566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&B)); 4549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N)); 4559566063dSJacob Faibussowitsch PetscCall(MatSetType(B,type)); 4569566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(B,ia,ja,a)); 45717667f90SBarry Smith 458511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 4599566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A,&B)); 46017667f90SBarry Smith } else { 46117667f90SBarry Smith *newmat = B; 46217667f90SBarry Smith } 46317667f90SBarry Smith PetscFunctionReturn(0); 46417667f90SBarry Smith } 46517667f90SBarry Smith 4666a09307cSBarry Smith PetscErrorCode MatSetValues_MPIAdj(Mat A,PetscInt m,const PetscInt *rows,PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode im) 4676a09307cSBarry Smith { 4686a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 4696a09307cSBarry Smith PetscInt rStart, rEnd, cStart, cEnd; 4706a09307cSBarry Smith 4716a09307cSBarry Smith PetscFunctionBegin; 4726a09307cSBarry Smith PetscCheck(!adj->i,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is already assembled, cannot change its values"); 4736a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 4746a09307cSBarry Smith PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd)); 4756a09307cSBarry Smith if (!adj->ht) { 4766a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht)); 4776a09307cSBarry Smith PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject) A), 1, &A->stash)); 4786a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 4796a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 4806a09307cSBarry Smith } 4816a09307cSBarry Smith for (PetscInt r = 0; r < m; ++r) { 4826a09307cSBarry Smith PetscHashIJKey key; 4836a09307cSBarry Smith 4846a09307cSBarry Smith key.i = rows[r]; 4856a09307cSBarry Smith if (key.i < 0) continue; 4866a09307cSBarry Smith if ((key.i < rStart) || (key.i >= rEnd)) { 4876a09307cSBarry Smith PetscCall(MatStashValuesRow_Private(&A->stash, key.i, n, cols, values, PETSC_FALSE)); 4886a09307cSBarry Smith } else { 4896a09307cSBarry Smith for (PetscInt c = 0; c < n; ++c) { 4906a09307cSBarry Smith key.j = cols[c]; 4916a09307cSBarry Smith if (key.j < 0 || key.i == key.j) continue; 4926a09307cSBarry Smith PetscCall(PetscHSetIJAdd(adj->ht, key)); 4936a09307cSBarry Smith } 4946a09307cSBarry Smith } 4956a09307cSBarry Smith } 4966a09307cSBarry Smith PetscFunctionReturn(0); 4976a09307cSBarry Smith } 4986a09307cSBarry Smith 4996a09307cSBarry Smith PetscErrorCode MatAssemblyBegin_MPIAdj(Mat A, MatAssemblyType type) 5006a09307cSBarry Smith { 5016a09307cSBarry Smith PetscInt nstash, reallocs; 5026a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 5036a09307cSBarry Smith 5046a09307cSBarry Smith PetscFunctionBegin; 5056a09307cSBarry Smith if (!adj->ht) { 5066a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht)); 5076a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 5086a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 5096a09307cSBarry Smith } 5106a09307cSBarry Smith PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range)); 5116a09307cSBarry Smith PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs)); 5126a09307cSBarry Smith PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs)); 5136a09307cSBarry Smith PetscFunctionReturn(0); 5146a09307cSBarry Smith } 5156a09307cSBarry Smith 5166a09307cSBarry Smith PetscErrorCode MatAssemblyEnd_MPIAdj(Mat A, MatAssemblyType type) 5176a09307cSBarry Smith { 5186a09307cSBarry Smith PetscScalar *val; 5196a09307cSBarry Smith PetscInt *row, *col,m,rstart,*rowstarts; 5206a09307cSBarry Smith PetscInt i, j, ncols, flg, nz; 5216a09307cSBarry Smith PetscMPIInt n; 5226a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 5236a09307cSBarry Smith PetscHashIter hi; 5246a09307cSBarry Smith PetscHashIJKey key; 5256a09307cSBarry Smith PetscHSetIJ ht = adj->ht; 5266a09307cSBarry Smith 5276a09307cSBarry Smith PetscFunctionBegin; 5286a09307cSBarry Smith while (1) { 5296a09307cSBarry Smith PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg)); 5306a09307cSBarry Smith if (!flg) break; 5316a09307cSBarry Smith 5326a09307cSBarry Smith for (i = 0; i < n;) { 5336a09307cSBarry Smith /* Identify the consecutive vals belonging to the same row */ 5346a09307cSBarry Smith for (j = i, rstart = row[j]; j < n; j++) { 5356a09307cSBarry Smith if (row[j] != rstart) break; 5366a09307cSBarry Smith } 5376a09307cSBarry Smith if (j < n) ncols = j-i; 5386a09307cSBarry Smith else ncols = n-i; 5396a09307cSBarry Smith /* Set all these values with a single function call */ 5406a09307cSBarry Smith PetscCall(MatSetValues_MPIAdj(A, 1, row+i, ncols, col+i, val+i, INSERT_VALUES)); 5416a09307cSBarry Smith i = j; 5426a09307cSBarry Smith } 5436a09307cSBarry Smith } 5446a09307cSBarry Smith PetscCall(MatStashScatterEnd_Private(&A->stash)); 5456a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&A->stash)); 5466a09307cSBarry Smith 5476a09307cSBarry Smith PetscCall(MatGetLocalSize(A,&m,NULL)); 5486a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A,&rstart,NULL)); 5496a09307cSBarry Smith PetscCall(PetscCalloc1(m+1,&rowstarts)); 5506a09307cSBarry Smith PetscHashIterBegin(ht,hi); 5516a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht,hi);) { 5526a09307cSBarry Smith PetscHashIterGetKey(ht,hi,key); 5536a09307cSBarry Smith rowstarts[key.i - rstart + 1]++; 5546a09307cSBarry Smith PetscHashIterNext(ht,hi); 5556a09307cSBarry Smith } 5566a09307cSBarry Smith for (i=1; i<m+1; i++) { 5576a09307cSBarry Smith rowstarts[i] = rowstarts[i-1] + rowstarts[i]; 5586a09307cSBarry Smith } 5596a09307cSBarry Smith 5606a09307cSBarry Smith PetscCall(PetscHSetIJGetSize(ht,&nz)); 5616a09307cSBarry Smith PetscCall(PetscMalloc1(nz,&col)); 5626a09307cSBarry Smith PetscHashIterBegin(ht,hi); 5636a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht,hi);) { 5646a09307cSBarry Smith PetscHashIterGetKey(ht,hi,key); 5656a09307cSBarry Smith col[rowstarts[key.i - rstart]++] = key.j; 5666a09307cSBarry Smith PetscHashIterNext(ht,hi); 5676a09307cSBarry Smith } 5686a09307cSBarry Smith PetscCall(PetscHSetIJDestroy(&ht)); 5696a09307cSBarry Smith for (i=m; i>0; i--) { 5706a09307cSBarry Smith rowstarts[i] = rowstarts[i-1]; 5716a09307cSBarry Smith } 5726a09307cSBarry Smith rowstarts[0] = 0; 5736a09307cSBarry Smith 5746a09307cSBarry Smith for (PetscInt i=0; i<m; i++) { 5756a09307cSBarry Smith PetscCall(PetscSortInt(rowstarts[i+1]-rowstarts[i],&col[rowstarts[i]])); 5766a09307cSBarry Smith } 5776a09307cSBarry Smith 5786a09307cSBarry Smith adj->i = rowstarts; 5796a09307cSBarry Smith adj->j = col; 5806a09307cSBarry Smith adj->freeaij = PETSC_TRUE; 5816a09307cSBarry Smith PetscFunctionReturn(0); 5826a09307cSBarry Smith } 5836a09307cSBarry Smith 584b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 5856a09307cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIAdj, 5863eda8832SBarry Smith MatGetRow_MPIAdj, 5873eda8832SBarry Smith MatRestoreRow_MPIAdj, 588f4259b30SLisandro Dalcin NULL, 589f4259b30SLisandro Dalcin /* 4*/ NULL, 590f4259b30SLisandro Dalcin NULL, 591f4259b30SLisandro Dalcin NULL, 592f4259b30SLisandro Dalcin NULL, 593f4259b30SLisandro Dalcin NULL, 594f4259b30SLisandro Dalcin NULL, 595f4259b30SLisandro Dalcin /*10*/ NULL, 596f4259b30SLisandro Dalcin NULL, 597f4259b30SLisandro Dalcin NULL, 598f4259b30SLisandro Dalcin NULL, 599f4259b30SLisandro Dalcin NULL, 600f4259b30SLisandro Dalcin /*15*/ NULL, 6013eda8832SBarry Smith MatEqual_MPIAdj, 602f4259b30SLisandro Dalcin NULL, 603f4259b30SLisandro Dalcin NULL, 604f4259b30SLisandro Dalcin NULL, 6056a09307cSBarry Smith /*20*/ MatAssemblyBegin_MPIAdj, 6066a09307cSBarry Smith MatAssemblyEnd_MPIAdj, 6073eda8832SBarry Smith MatSetOption_MPIAdj, 608f4259b30SLisandro Dalcin NULL, 609f4259b30SLisandro Dalcin /*24*/ NULL, 610f4259b30SLisandro Dalcin NULL, 611f4259b30SLisandro Dalcin NULL, 612f4259b30SLisandro Dalcin NULL, 613f4259b30SLisandro Dalcin NULL, 614f4259b30SLisandro Dalcin /*29*/ NULL, 615f4259b30SLisandro Dalcin NULL, 616f4259b30SLisandro Dalcin NULL, 617f4259b30SLisandro Dalcin NULL, 618f4259b30SLisandro Dalcin NULL, 619f4259b30SLisandro Dalcin /*34*/ NULL, 620f4259b30SLisandro Dalcin NULL, 621f4259b30SLisandro Dalcin NULL, 622f4259b30SLisandro Dalcin NULL, 623f4259b30SLisandro Dalcin NULL, 624f4259b30SLisandro Dalcin /*39*/ NULL, 6257dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 626f4259b30SLisandro Dalcin NULL, 627f4259b30SLisandro Dalcin NULL, 628f4259b30SLisandro Dalcin NULL, 629f4259b30SLisandro Dalcin /*44*/ NULL, 630f4259b30SLisandro Dalcin NULL, 6317d68702bSBarry Smith MatShift_Basic, 632f4259b30SLisandro Dalcin NULL, 633f4259b30SLisandro Dalcin NULL, 634f4259b30SLisandro Dalcin /*49*/ NULL, 6356cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 636d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 637f4259b30SLisandro Dalcin NULL, 638f4259b30SLisandro Dalcin NULL, 639f4259b30SLisandro Dalcin /*54*/ NULL, 640f4259b30SLisandro Dalcin NULL, 641f4259b30SLisandro Dalcin NULL, 642f4259b30SLisandro Dalcin NULL, 643f4259b30SLisandro Dalcin NULL, 644f4259b30SLisandro Dalcin /*59*/ NULL, 645b9b97703SBarry Smith MatDestroy_MPIAdj, 646b9b97703SBarry Smith MatView_MPIAdj, 64717667f90SBarry Smith MatConvertFrom_MPIAdj, 648f4259b30SLisandro Dalcin NULL, 649f4259b30SLisandro Dalcin /*64*/ NULL, 650f4259b30SLisandro Dalcin NULL, 651f4259b30SLisandro Dalcin NULL, 652f4259b30SLisandro Dalcin NULL, 653f4259b30SLisandro Dalcin NULL, 654f4259b30SLisandro Dalcin /*69*/ NULL, 655f4259b30SLisandro Dalcin NULL, 656f4259b30SLisandro Dalcin NULL, 657f4259b30SLisandro Dalcin NULL, 658f4259b30SLisandro Dalcin NULL, 659f4259b30SLisandro Dalcin /*74*/ NULL, 660f4259b30SLisandro Dalcin NULL, 661f4259b30SLisandro Dalcin NULL, 662f4259b30SLisandro Dalcin NULL, 663f4259b30SLisandro Dalcin NULL, 664f4259b30SLisandro Dalcin /*79*/ NULL, 665f4259b30SLisandro Dalcin NULL, 666f4259b30SLisandro Dalcin NULL, 667f4259b30SLisandro Dalcin NULL, 668f4259b30SLisandro Dalcin NULL, 669f4259b30SLisandro Dalcin /*84*/ NULL, 670f4259b30SLisandro Dalcin NULL, 671f4259b30SLisandro Dalcin NULL, 672f4259b30SLisandro Dalcin NULL, 673f4259b30SLisandro Dalcin NULL, 674f4259b30SLisandro Dalcin /*89*/ NULL, 675f4259b30SLisandro Dalcin NULL, 676f4259b30SLisandro Dalcin NULL, 677f4259b30SLisandro Dalcin NULL, 678f4259b30SLisandro Dalcin NULL, 679f4259b30SLisandro Dalcin /*94*/ NULL, 680f4259b30SLisandro Dalcin NULL, 681f4259b30SLisandro Dalcin NULL, 682f4259b30SLisandro Dalcin NULL, 683f4259b30SLisandro Dalcin NULL, 684f4259b30SLisandro Dalcin /*99*/ NULL, 685f4259b30SLisandro Dalcin NULL, 686f4259b30SLisandro Dalcin NULL, 687f4259b30SLisandro Dalcin NULL, 688f4259b30SLisandro Dalcin NULL, 689f4259b30SLisandro Dalcin /*104*/ NULL, 690f4259b30SLisandro Dalcin NULL, 691f4259b30SLisandro Dalcin NULL, 692f4259b30SLisandro Dalcin NULL, 693f4259b30SLisandro Dalcin NULL, 694f4259b30SLisandro Dalcin /*109*/ NULL, 695f4259b30SLisandro Dalcin NULL, 696f4259b30SLisandro Dalcin NULL, 697f4259b30SLisandro Dalcin NULL, 698f4259b30SLisandro Dalcin NULL, 699f4259b30SLisandro Dalcin /*114*/ NULL, 700f4259b30SLisandro Dalcin NULL, 701f4259b30SLisandro Dalcin NULL, 702f4259b30SLisandro Dalcin NULL, 703f4259b30SLisandro Dalcin NULL, 704f4259b30SLisandro Dalcin /*119*/ NULL, 705f4259b30SLisandro Dalcin NULL, 706f4259b30SLisandro Dalcin NULL, 707f4259b30SLisandro Dalcin NULL, 708f4259b30SLisandro Dalcin NULL, 709f4259b30SLisandro Dalcin /*124*/ NULL, 710f4259b30SLisandro Dalcin NULL, 711f4259b30SLisandro Dalcin NULL, 712f4259b30SLisandro Dalcin NULL, 7137dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 714f4259b30SLisandro Dalcin /*129*/ NULL, 715f4259b30SLisandro Dalcin NULL, 716f4259b30SLisandro Dalcin NULL, 717f4259b30SLisandro Dalcin NULL, 718f4259b30SLisandro Dalcin NULL, 719f4259b30SLisandro Dalcin /*134*/ NULL, 720f4259b30SLisandro Dalcin NULL, 721f4259b30SLisandro Dalcin NULL, 722f4259b30SLisandro Dalcin NULL, 723f4259b30SLisandro Dalcin NULL, 724f4259b30SLisandro Dalcin /*139*/ NULL, 725f4259b30SLisandro Dalcin NULL, 726d70f29a3SPierre Jolivet NULL, 727d70f29a3SPierre Jolivet NULL, 728d70f29a3SPierre Jolivet NULL, 729d70f29a3SPierre Jolivet /*144*/NULL, 730d70f29a3SPierre Jolivet NULL, 731d70f29a3SPierre Jolivet NULL, 73299a7f59eSMark Adams NULL, 73399a7f59eSMark Adams NULL, 734f4259b30SLisandro Dalcin NULL 7353964eb88SJed Brown }; 736b97cf49bSBarry Smith 737f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 738a23d5eceSKris Buschelman { 739a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 740e895ccc0SFande Kong PetscBool useedgeweights; 741a23d5eceSKris Buschelman 742a23d5eceSKris Buschelman PetscFunctionBegin; 7439566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->rmap)); 7449566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->cmap)); 745e895ccc0SFande Kong if (values) useedgeweights = PETSC_TRUE; else useedgeweights = PETSC_FALSE; 746e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */ 7471c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((int*)&useedgeweights,(int*)&b->useedgeweights,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)B))); 74858ec18a5SLisandro Dalcin 74976bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 75076bd3646SJed Brown PetscInt ii; 75176bd3646SJed Brown 75208401ef6SPierre 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]); 753d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 754aed4548fSBarry 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]); 755a23d5eceSKris Buschelman } 756d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 757aed4548fSBarry 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]); 758a23d5eceSKris Buschelman } 75976bd3646SJed Brown } 760a23d5eceSKris Buschelman b->j = j; 761a23d5eceSKris Buschelman b->i = i; 762a23d5eceSKris Buschelman b->values = values; 763a23d5eceSKris Buschelman 764d0f46423SBarry Smith b->nz = i[B->rmap->n]; 765f4259b30SLisandro Dalcin b->diag = NULL; 766a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 767a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 768a23d5eceSKris Buschelman 7696a09307cSBarry Smith B->ops->assemblybegin = NULL; 7706a09307cSBarry Smith B->ops->assemblyend = NULL; 7719566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 7729566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 7736a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&B->stash)); 774a23d5eceSKris Buschelman PetscFunctionReturn(0); 775a23d5eceSKris Buschelman } 776b97cf49bSBarry Smith 777f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 7789a3665c6SJed Brown { 7799a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 7809a3665c6SJed Brown const PetscInt *ranges; 7819a3665c6SJed Brown MPI_Comm acomm,bcomm; 7829a3665c6SJed Brown MPI_Group agroup,bgroup; 7839a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 7849a3665c6SJed Brown 7859a3665c6SJed Brown PetscFunctionBegin; 7860298fd71SBarry Smith *B = NULL; 7879566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A,&acomm)); 7889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm,&size)); 7899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm,&rank)); 7909566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A,&ranges)); 7919a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 7929a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 7939a3665c6SJed Brown } 7949a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 7959566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 7969a3665c6SJed Brown *B = A; 7979a3665c6SJed Brown PetscFunctionReturn(0); 7989a3665c6SJed Brown } 7999a3665c6SJed Brown 8009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nranks,&ranks)); 8019a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 8029a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 8039a3665c6SJed Brown } 8049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(acomm,&agroup)); 8059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(agroup,nranks,ranks,&bgroup)); 8069566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks)); 8079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(acomm,bgroup,&bcomm)); 8089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&agroup)); 8099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&bgroup)); 8109a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 8119a3665c6SJed Brown PetscInt m,N; 8129a3665c6SJed Brown Mat_MPIAdj *b; 8139566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 8149566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 8159566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B)); 8169a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 8179a3665c6SJed Brown b->freeaij = PETSC_FALSE; 8189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&bcomm)); 8199a3665c6SJed Brown } 8209a3665c6SJed Brown PetscFunctionReturn(0); 8219a3665c6SJed Brown } 8229a3665c6SJed Brown 823b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B) 824b44e7856SBarry Smith { 825b44e7856SBarry Smith PetscInt M,N,*II,*J,NZ,nz,m,nzstart,i; 826b44e7856SBarry Smith PetscInt *Values = NULL; 827b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 828b44e7856SBarry Smith PetscMPIInt mnz,mm,*allnz,*allm,size,*dispnz,*dispm; 829b44e7856SBarry Smith 830b44e7856SBarry Smith PetscFunctionBegin; 8319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 8329566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 8339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 834b44e7856SBarry Smith nz = adj->nz; 83508401ef6SPierre 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]); 8369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A))); 837d4e69552SBarry Smith 8389566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nz,&mnz)); 8399566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size,&allnz,size,&dispnz)); 8409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A))); 841b44e7856SBarry Smith dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1]; 842b44e7856SBarry Smith if (adj->values) { 8439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ,&Values)); 8449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A))); 845b44e7856SBarry Smith } 8469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ,&J)); 8479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A))); 8489566063dSJacob Faibussowitsch PetscCall(PetscFree2(allnz,dispnz)); 8499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A))); 850b44e7856SBarry Smith nzstart -= nz; 851b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 852b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] += nzstart; 8539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M+1,&II)); 8549566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(m,&mm)); 8559566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size,&allm,size,&dispm)); 8569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A))); 857b44e7856SBarry Smith dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1]; 8589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A))); 8599566063dSJacob Faibussowitsch PetscCall(PetscFree2(allm,dispm)); 860b44e7856SBarry Smith II[M] = NZ; 861b44e7856SBarry Smith /* shift the i[] values back */ 862b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] -= nzstart; 8639566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B)); 864b44e7856SBarry Smith PetscFunctionReturn(0); 865b44e7856SBarry Smith } 866b44e7856SBarry Smith 8679a3665c6SJed Brown /*@ 8689a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 8699a3665c6SJed Brown 8709a3665c6SJed Brown Collective 8719a3665c6SJed Brown 8724165533cSJose E. Roman Input Parameter: 8739a3665c6SJed Brown . A - original MPIAdj matrix 8749a3665c6SJed Brown 8754165533cSJose E. Roman Output Parameter: 8760298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 8779a3665c6SJed Brown 8789a3665c6SJed Brown Level: developer 8799a3665c6SJed Brown 8809a3665c6SJed Brown Note: 8819a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 8829a3665c6SJed Brown 8839a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 8849a3665c6SJed Brown 885db781477SPatrick Sanan .seealso: `MatCreateMPIAdj()` 8869a3665c6SJed Brown @*/ 8879a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 8889a3665c6SJed Brown { 8899a3665c6SJed Brown PetscFunctionBegin; 8909a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 891cac4c232SBarry Smith PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B)); 8929a3665c6SJed Brown PetscFunctionReturn(0); 8939a3665c6SJed Brown } 8949a3665c6SJed Brown 8950bad9183SKris Buschelman /*MC 896fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 8970bad9183SKris Buschelman intended for use constructing orderings and partitionings. 8980bad9183SKris Buschelman 8990bad9183SKris Buschelman Level: beginner 9000bad9183SKris Buschelman 9016a09307cSBarry Smith Notes: 9026a09307cSBarry Smith You can provide values to the matrix using MatMPIAdjSetPreallocation(), MatCreateMPIAdj(), or 9036a09307cSBarry Smith by calling MatSetValues() and MatAssemblyBegin() followed by MatAssemblyEnd() 9046a09307cSBarry Smith 9056a09307cSBarry Smith .seealso: `MatCreateMPIAdj()`, `MatMPIAdjSetPreallocation()`, `MatSetValues()` 9060bad9183SKris Buschelman M*/ 9070bad9183SKris Buschelman 9088cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 909273d9f13SBarry Smith { 910273d9f13SBarry Smith Mat_MPIAdj *b; 911273d9f13SBarry Smith 912273d9f13SBarry Smith PetscFunctionBegin; 9139566063dSJacob Faibussowitsch PetscCall(PetscNewLog(B,&b)); 914b0a32e0cSBarry Smith B->data = (void*)b; 9159566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps))); 916273d9f13SBarry Smith B->assembled = PETSC_FALSE; 9176a09307cSBarry Smith B->preallocated = PETSC_TRUE; /* so that MatSetValues() may be used */ 918273d9f13SBarry Smith 9199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj)); 9209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj)); 9219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj)); 9229566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ)); 923273d9f13SBarry Smith PetscFunctionReturn(0); 924273d9f13SBarry Smith } 925273d9f13SBarry Smith 926a23d5eceSKris Buschelman /*@C 927b44e7856SBarry Smith MatMPIAdjToSeq - Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential preconditioners) 928b44e7856SBarry Smith 929d083f849SBarry Smith Logically Collective 930b44e7856SBarry Smith 931b44e7856SBarry Smith Input Parameter: 932b44e7856SBarry Smith . A - the matrix 933b44e7856SBarry Smith 934b44e7856SBarry Smith Output Parameter: 935b44e7856SBarry Smith . B - the same matrix on all processes 936b44e7856SBarry Smith 937b44e7856SBarry Smith Level: intermediate 938b44e7856SBarry Smith 939db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()` 940b44e7856SBarry Smith @*/ 941b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq(Mat A,Mat *B) 942b44e7856SBarry Smith { 943b44e7856SBarry Smith PetscFunctionBegin; 944cac4c232SBarry Smith PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B)); 945b44e7856SBarry Smith PetscFunctionReturn(0); 946b44e7856SBarry Smith } 947b44e7856SBarry Smith 948b44e7856SBarry Smith /*@C 949a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 950a23d5eceSKris Buschelman 951d083f849SBarry Smith Logically Collective 952a23d5eceSKris Buschelman 953a23d5eceSKris Buschelman Input Parameters: 954a23d5eceSKris Buschelman + A - the matrix 955a23d5eceSKris Buschelman . i - the indices into j for the start of each row 956a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 957a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 958a23d5eceSKris Buschelman - values - [optional] edge weights 959a23d5eceSKris Buschelman 960a23d5eceSKris Buschelman Level: intermediate 961a23d5eceSKris Buschelman 9626a09307cSBarry Smith .seealso: `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MATMPIADJ`, `MatCreateMPIAdj()` 963a23d5eceSKris Buschelman @*/ 9647087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 965273d9f13SBarry Smith { 966273d9f13SBarry Smith PetscFunctionBegin; 967cac4c232SBarry Smith PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values)); 968273d9f13SBarry Smith PetscFunctionReturn(0); 969273d9f13SBarry Smith } 970273d9f13SBarry Smith 971b97cf49bSBarry Smith /*@C 9723eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 973b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 974b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 975b97cf49bSBarry Smith 976d083f849SBarry Smith Collective 977ef5ee4d1SLois Curfman McInnes 978b97cf49bSBarry Smith Input Parameters: 979c2e958feSBarry Smith + comm - MPI communicator 9800752156aSBarry Smith . m - number of local rows 98158ec18a5SLisandro Dalcin . N - number of global columns 982b97cf49bSBarry Smith . i - the indices into j for the start of each row 983329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 984ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 985329f5518SBarry Smith - values -[optional] edge weights 986b97cf49bSBarry Smith 987b97cf49bSBarry Smith Output Parameter: 988b97cf49bSBarry Smith . A - the matrix 989b97cf49bSBarry Smith 99015091d37SBarry Smith Level: intermediate 99115091d37SBarry Smith 99295452b02SPatrick Sanan Notes: 993c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 994ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 9951198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 9966a09307cSBarry Smith 9976a09307cSBarry Smith You should not include the matrix diagonals. 998b97cf49bSBarry Smith 999e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 10006a09307cSBarry Smith to MatConvert(), specifying a type of MATMPIADJ. 1001ededeb1aSvictorle 1002ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 1003b97cf49bSBarry Smith 10046a09307cSBarry Smith .seealso: `MatCreate()`, `MatConvert()`, `MatGetOrdering()`, `MATMPIADJ`, `MatMPIAdjSetPreallocation()` 1005b97cf49bSBarry Smith @*/ 10067087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 1007b97cf49bSBarry Smith { 1008433994e6SBarry Smith PetscFunctionBegin; 10099566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,A)); 10109566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N)); 10119566063dSJacob Faibussowitsch PetscCall(MatSetType(*A,MATMPIADJ)); 10129566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(*A,i,j,values)); 10133a40ed3dSBarry Smith PetscFunctionReturn(0); 1014b97cf49bSBarry Smith } 1015