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 { 13*131c27b5Sprj- 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; 17*131c27b5Sprj- 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 PetscErrorCode ierr; 2540244768SBarry Smith 2640244768SBarry Smith PetscFunctionBegin; 2740244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)adj,&comm);CHKERRQ(ierr); 28390e1bf2SBarry Smith ierr = MatGetLayouts(adj,&rmap,NULL);CHKERRQ(ierr); 2940244768SBarry Smith ierr = ISGetLocalSize(irows,&nlrows_is);CHKERRQ(ierr); 3040244768SBarry Smith ierr = ISGetIndices(irows,&irows_indices);CHKERRQ(ierr); 31071fcb05SBarry Smith ierr = PetscMalloc1(nlrows_is,&iremote);CHKERRQ(ierr); 3240244768SBarry Smith /* construct sf graph*/ 3340244768SBarry Smith nleaves = nlrows_is; 3440244768SBarry Smith for (i=0; i<nlrows_is; i++){ 3540244768SBarry Smith owner = -1; 3640244768SBarry Smith rlocalindex = -1; 3740244768SBarry Smith ierr = PetscLayoutFindOwnerIndex(rmap,irows_indices[i],&owner,&rlocalindex);CHKERRQ(ierr); 3840244768SBarry Smith iremote[i].rank = owner; 3940244768SBarry Smith iremote[i].index = rlocalindex; 4040244768SBarry Smith } 4140244768SBarry Smith ierr = MatGetRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr); 4240244768SBarry Smith ierr = PetscCalloc4(nlrows_mat,&ncols_send,nlrows_is,&xadj_recv,nlrows_is+1,&ncols_recv_offsets,nlrows_is,&ncols_recv);CHKERRQ(ierr); 4340244768SBarry Smith nroots = nlrows_mat; 4440244768SBarry Smith for (i=0; i<nlrows_mat; i++){ 4540244768SBarry Smith ncols_send[i] = xadj[i+1]-xadj[i]; 4640244768SBarry Smith } 4740244768SBarry Smith ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 48390e1bf2SBarry Smith ierr = PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 4940244768SBarry Smith ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 5040244768SBarry Smith ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 5140244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr); 5240244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr); 5340244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr); 5440244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr); 5540244768SBarry Smith ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 5640244768SBarry Smith Ncols_recv =0; 5740244768SBarry Smith for (i=0; i<nlrows_is; i++){ 5840244768SBarry Smith Ncols_recv += ncols_recv[i]; 5940244768SBarry Smith ncols_recv_offsets[i+1] = ncols_recv[i]+ncols_recv_offsets[i]; 6040244768SBarry Smith } 6140244768SBarry Smith Ncols_send = 0; 6240244768SBarry Smith for(i=0; i<nlrows_mat; i++){ 6340244768SBarry Smith Ncols_send += ncols_send[i]; 6440244768SBarry Smith } 6540244768SBarry Smith ierr = PetscCalloc1(Ncols_recv,&iremote);CHKERRQ(ierr); 6640244768SBarry Smith ierr = PetscCalloc1(Ncols_recv,&adjncy_recv);CHKERRQ(ierr); 6740244768SBarry Smith nleaves = Ncols_recv; 6840244768SBarry Smith Ncols_recv = 0; 6940244768SBarry Smith for (i=0; i<nlrows_is; i++){ 7040244768SBarry Smith ierr = PetscLayoutFindOwner(rmap,irows_indices[i],&owner);CHKERRQ(ierr); 7140244768SBarry Smith for (j=0; j<ncols_recv[i]; j++){ 7240244768SBarry Smith iremote[Ncols_recv].rank = owner; 7340244768SBarry Smith iremote[Ncols_recv++].index = xadj_recv[i]+j; 7440244768SBarry Smith } 7540244768SBarry Smith } 7640244768SBarry Smith ierr = ISRestoreIndices(irows,&irows_indices);CHKERRQ(ierr); 7740244768SBarry Smith /*if we need to deal with edge weights ???*/ 78e895ccc0SFande Kong if (a->useedgeweights){ierr = PetscCalloc1(Ncols_recv,&values_recv);CHKERRQ(ierr);} 7940244768SBarry Smith nroots = Ncols_send; 8040244768SBarry Smith ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 81390e1bf2SBarry Smith ierr = PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 8240244768SBarry Smith ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 8340244768SBarry Smith ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 8440244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 8540244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 86e895ccc0SFande Kong if (a->useedgeweights){ 8740244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 8840244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 8940244768SBarry Smith } 9040244768SBarry Smith ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 9140244768SBarry Smith ierr = MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr); 9240244768SBarry Smith ierr = ISGetLocalSize(icols,&icols_n);CHKERRQ(ierr); 9340244768SBarry Smith ierr = ISGetIndices(icols,&icols_indices);CHKERRQ(ierr); 9440244768SBarry Smith rnclos = 0; 9540244768SBarry Smith for (i=0; i<nlrows_is; i++){ 9640244768SBarry Smith for (j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 9740244768SBarry Smith ierr = PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc);CHKERRQ(ierr); 9840244768SBarry Smith if (loc<0){ 9940244768SBarry Smith adjncy_recv[j] = -1; 100e895ccc0SFande Kong if (a->useedgeweights) values_recv[j] = -1; 10140244768SBarry Smith ncols_recv[i]--; 10240244768SBarry Smith } else { 10340244768SBarry Smith rnclos++; 10440244768SBarry Smith } 10540244768SBarry Smith } 10640244768SBarry Smith } 10740244768SBarry Smith ierr = ISRestoreIndices(icols,&icols_indices);CHKERRQ(ierr); 10840244768SBarry Smith ierr = PetscCalloc1(rnclos,&sadjncy);CHKERRQ(ierr); 109e895ccc0SFande Kong if (a->useedgeweights) {ierr = PetscCalloc1(rnclos,&svalues);CHKERRQ(ierr);} 11040244768SBarry Smith ierr = PetscCalloc1(nlrows_is+1,&sxadj);CHKERRQ(ierr); 11140244768SBarry Smith rnclos = 0; 11240244768SBarry Smith for(i=0; i<nlrows_is; i++){ 11340244768SBarry Smith for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 11440244768SBarry Smith if (adjncy_recv[j]<0) continue; 11540244768SBarry Smith sadjncy[rnclos] = adjncy_recv[j]; 116e895ccc0SFande Kong if (a->useedgeweights) svalues[rnclos] = values_recv[j]; 11740244768SBarry Smith rnclos++; 11840244768SBarry Smith } 11940244768SBarry Smith } 12040244768SBarry Smith for(i=0; i<nlrows_is; i++){ 12140244768SBarry Smith sxadj[i+1] = sxadj[i]+ncols_recv[i]; 12240244768SBarry Smith } 12340244768SBarry Smith if (sadj_xadj) { *sadj_xadj = sxadj;} else { ierr = PetscFree(sxadj);CHKERRQ(ierr);} 12440244768SBarry Smith if (sadj_adjncy){ *sadj_adjncy = sadjncy;} else { ierr = PetscFree(sadjncy);CHKERRQ(ierr);} 12540244768SBarry Smith if (sadj_values){ 126e895ccc0SFande Kong if (a->useedgeweights) *sadj_values = svalues; else *sadj_values=0; 12740244768SBarry Smith } else { 128e895ccc0SFande Kong if (a->useedgeweights) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 12940244768SBarry Smith } 13040244768SBarry Smith ierr = PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv);CHKERRQ(ierr); 13140244768SBarry Smith ierr = PetscFree(adjncy_recv);CHKERRQ(ierr); 132e895ccc0SFande Kong if (a->useedgeweights) {ierr = PetscFree(values_recv);CHKERRQ(ierr);} 13340244768SBarry Smith PetscFunctionReturn(0); 13440244768SBarry Smith } 13540244768SBarry Smith 1367dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj_Private(Mat mat,PetscInt n,const IS irow[],const IS icol[],PetscBool subcomm,MatReuse scall,Mat *submat[]) 13740244768SBarry Smith { 13840244768SBarry Smith PetscInt i,irow_n,icol_n,*sxadj,*sadjncy,*svalues; 13940244768SBarry Smith PetscInt *indices,nindx,j,k,loc; 14040244768SBarry Smith PetscMPIInt issame; 14140244768SBarry Smith const PetscInt *irow_indices,*icol_indices; 14240244768SBarry Smith MPI_Comm scomm_row,scomm_col,scomm_mat; 14340244768SBarry Smith PetscErrorCode ierr; 14440244768SBarry Smith 14540244768SBarry Smith PetscFunctionBegin; 14640244768SBarry Smith nindx = 0; 14740244768SBarry Smith /* 14840244768SBarry Smith * Estimate a maximum number for allocating memory 14940244768SBarry Smith */ 15040244768SBarry Smith for(i=0; i<n; i++){ 15140244768SBarry Smith ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr); 15240244768SBarry Smith ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr); 15340244768SBarry Smith nindx = nindx>(irow_n+icol_n)? nindx:(irow_n+icol_n); 15440244768SBarry Smith } 155071fcb05SBarry Smith ierr = PetscMalloc1(nindx,&indices);CHKERRQ(ierr); 15640244768SBarry Smith /* construct a submat */ 15740244768SBarry Smith for (i=0; i<n; i++){ 15840244768SBarry Smith if (subcomm){ 15940244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)irow[i],&scomm_row);CHKERRQ(ierr); 16040244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)icol[i],&scomm_col);CHKERRQ(ierr); 16140244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,scomm_col,&issame);CHKERRQ(ierr); 16240244768SBarry Smith if (issame != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row index set must have the same comm as the col index set\n"); 16340244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame);CHKERRQ(ierr); 16440244768SBarry Smith if (issame == MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP," can not use PETSC_COMM_SELF as comm when extracting a parallel submatrix\n"); 16540244768SBarry Smith } else { 16640244768SBarry Smith scomm_row = PETSC_COMM_SELF; 16740244768SBarry Smith } 16840244768SBarry Smith /*get sub-matrix data*/ 16940244768SBarry Smith sxadj=0; sadjncy=0; svalues=0; 1707dae84e0SHong Zhang ierr = MatCreateSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues);CHKERRQ(ierr); 17140244768SBarry Smith ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr); 17240244768SBarry Smith ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr); 17340244768SBarry Smith ierr = ISGetIndices(irow[i],&irow_indices);CHKERRQ(ierr); 174580bdb30SBarry Smith ierr = PetscArraycpy(indices,irow_indices,irow_n);CHKERRQ(ierr); 17540244768SBarry Smith ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr); 17640244768SBarry Smith ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr); 177580bdb30SBarry Smith ierr = PetscArraycpy(indices+irow_n,icol_indices,icol_n);CHKERRQ(ierr); 17840244768SBarry Smith ierr = ISRestoreIndices(icol[i],&icol_indices);CHKERRQ(ierr); 17940244768SBarry Smith nindx = irow_n+icol_n; 18040244768SBarry Smith ierr = PetscSortRemoveDupsInt(&nindx,indices);CHKERRQ(ierr); 18140244768SBarry Smith /* renumber columns */ 18240244768SBarry Smith for (j=0; j<irow_n; j++){ 18340244768SBarry Smith for (k=sxadj[j]; k<sxadj[j+1]; k++){ 18440244768SBarry Smith ierr = PetscFindInt(sadjncy[k],nindx,indices,&loc);CHKERRQ(ierr); 18562795ce3SStefano Zampini #if defined(PETSC_USE_DEBUG) 18662795ce3SStefano Zampini if (loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %D",sadjncy[k]); 18740244768SBarry Smith #endif 18840244768SBarry Smith sadjncy[k] = loc; 18940244768SBarry Smith } 19040244768SBarry Smith } 19140244768SBarry Smith if (scall==MAT_INITIAL_MATRIX){ 19240244768SBarry Smith ierr = MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr); 19340244768SBarry Smith } else { 19440244768SBarry Smith Mat sadj = *(submat[i]); 19540244768SBarry Smith Mat_MPIAdj *sa = (Mat_MPIAdj*)((sadj)->data); 19640244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)sadj,&scomm_mat);CHKERRQ(ierr); 19740244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,scomm_mat,&issame);CHKERRQ(ierr); 19840244768SBarry Smith if (issame != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"submatrix must have the same comm as the col index set\n"); 199580bdb30SBarry Smith ierr = PetscArraycpy(sa->i,sxadj,irow_n+1);CHKERRQ(ierr); 200580bdb30SBarry Smith ierr = PetscArraycpy(sa->j,sadjncy,sxadj[irow_n]);CHKERRQ(ierr); 201580bdb30SBarry Smith if (svalues){ierr = PetscArraycpy(sa->values,svalues,sxadj[irow_n]);CHKERRQ(ierr);} 20240244768SBarry Smith ierr = PetscFree(sxadj);CHKERRQ(ierr); 20340244768SBarry Smith ierr = PetscFree(sadjncy);CHKERRQ(ierr); 20440244768SBarry Smith if (svalues) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 20540244768SBarry Smith } 20640244768SBarry Smith } 20740244768SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 20840244768SBarry Smith PetscFunctionReturn(0); 20940244768SBarry Smith } 21040244768SBarry Smith 2117dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatricesMPI_MPIAdj(Mat mat,PetscInt n, const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 21240244768SBarry Smith { 21340244768SBarry Smith PetscErrorCode ierr; 21440244768SBarry Smith /*get sub-matrices across a sub communicator */ 21540244768SBarry Smith PetscFunctionBegin; 2167dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat);CHKERRQ(ierr); 21740244768SBarry Smith PetscFunctionReturn(0); 21840244768SBarry Smith } 21940244768SBarry Smith 2207dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 22140244768SBarry Smith { 22240244768SBarry Smith PetscErrorCode ierr; 22340244768SBarry Smith 22440244768SBarry Smith PetscFunctionBegin; 22540244768SBarry Smith /*get sub-matrices based on PETSC_COMM_SELF */ 2267dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat);CHKERRQ(ierr); 22740244768SBarry Smith PetscFunctionReturn(0); 22840244768SBarry Smith } 22940244768SBarry Smith 23040244768SBarry Smith static PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 231b97cf49bSBarry Smith { 2323eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 233dfbe8321SBarry Smith PetscErrorCode ierr; 234d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 2352dcb1b2aSMatthew Knepley const char *name; 236ce1411ecSBarry Smith PetscViewerFormat format; 237b97cf49bSBarry Smith 238433994e6SBarry Smith PetscFunctionBegin; 2393a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 240b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 241fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 2423a40ed3dSBarry Smith PetscFunctionReturn(0); 2436c4ed002SBarry Smith } else if (format == PETSC_VIEWER_ASCII_MATLAB) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported"); 2446c4ed002SBarry Smith else { 245d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 246c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 247b97cf49bSBarry Smith for (i=0; i<m; i++) { 248d0f46423SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr); 249b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 2500170919eSFande Kong if (a->values) { 2510170919eSFande Kong ierr = PetscViewerASCIISynchronizedPrintf(viewer," (%D, %D) ",a->j[j], a->values[j]);CHKERRQ(ierr); 2520170919eSFande Kong } else { 25377431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 2540752156aSBarry Smith } 2550170919eSFande Kong } 256b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 257b97cf49bSBarry Smith } 258d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 259b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 260c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2617b23a99aSBarry Smith } 2623a40ed3dSBarry Smith PetscFunctionReturn(0); 263b97cf49bSBarry Smith } 264b97cf49bSBarry Smith 26540244768SBarry Smith static PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 266b97cf49bSBarry Smith { 267dfbe8321SBarry Smith PetscErrorCode ierr; 268ace3abfcSBarry Smith PetscBool iascii; 269b97cf49bSBarry Smith 270433994e6SBarry Smith PetscFunctionBegin; 271251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 27232077d6dSBarry Smith if (iascii) { 2733eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 274b97cf49bSBarry Smith } 2753a40ed3dSBarry Smith PetscFunctionReturn(0); 276b97cf49bSBarry Smith } 277b97cf49bSBarry Smith 27840244768SBarry Smith static PetscErrorCode MatDestroy_MPIAdj(Mat mat) 279b97cf49bSBarry Smith { 2803eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 281dfbe8321SBarry Smith PetscErrorCode ierr; 28294d884c6SBarry Smith 28394d884c6SBarry Smith PetscFunctionBegin; 284aa482453SBarry Smith #if defined(PETSC_USE_LOG) 285d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz); 286b97cf49bSBarry Smith #endif 28705b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 2880400557aSBarry Smith if (a->freeaij) { 28914196391SBarry Smith if (a->freeaijwithfree) { 29014196391SBarry Smith if (a->i) free(a->i); 29114196391SBarry Smith if (a->j) free(a->j); 29214196391SBarry Smith } else { 293606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 294606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 29505b42c5fSBarry Smith ierr = PetscFree(a->values);CHKERRQ(ierr); 29614196391SBarry Smith } 2970400557aSBarry Smith } 2982b1d8763SJed Brown ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 299bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 300dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 301bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr); 302bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr); 3033a40ed3dSBarry Smith PetscFunctionReturn(0); 304b97cf49bSBarry Smith } 305b97cf49bSBarry Smith 30640244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 307b97cf49bSBarry Smith { 3083eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 30963ba0a88SBarry Smith PetscErrorCode ierr; 310b97cf49bSBarry Smith 311433994e6SBarry Smith PetscFunctionBegin; 31212c028f9SKris Buschelman switch (op) { 31377e54ba9SKris Buschelman case MAT_SYMMETRIC: 31412c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3159a4540c5SBarry Smith case MAT_HERMITIAN: 3164e0d8c25SBarry Smith a->symmetric = flg; 3179a4540c5SBarry Smith break; 3189a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 3199a4540c5SBarry Smith break; 32012c028f9SKris Buschelman default: 321290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 32212c028f9SKris Buschelman break; 323b97cf49bSBarry Smith } 3243a40ed3dSBarry Smith PetscFunctionReturn(0); 325b97cf49bSBarry Smith } 326b97cf49bSBarry Smith 32740244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 328b97cf49bSBarry Smith { 3293eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3302b1d8763SJed Brown PetscErrorCode ierr; 331b97cf49bSBarry Smith 332433994e6SBarry Smith PetscFunctionBegin; 333d0f46423SBarry Smith row -= A->rmap->rstart; 334e32f2f54SBarry Smith if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 335b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 3362b1d8763SJed Brown if (v) { 3372b1d8763SJed Brown PetscInt j; 3382b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 3392b1d8763SJed Brown ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 3402b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz); 341785e854fSJed Brown ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr); 342b97cf49bSBarry Smith } 34391f8cafeSFande Kong for (j=0; j<*nz; j++){ 34491f8cafeSFande Kong a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0; 34591f8cafeSFande Kong } 3462b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 347b97cf49bSBarry Smith } 34892b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 3493a40ed3dSBarry Smith PetscFunctionReturn(0); 350b97cf49bSBarry Smith } 351b97cf49bSBarry Smith 35240244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 353b97cf49bSBarry Smith { 354433994e6SBarry Smith PetscFunctionBegin; 3553a40ed3dSBarry Smith PetscFunctionReturn(0); 356b97cf49bSBarry Smith } 357b97cf49bSBarry Smith 35840244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 359b97cf49bSBarry Smith { 3603eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 361dfbe8321SBarry Smith PetscErrorCode ierr; 362ace3abfcSBarry Smith PetscBool flag; 363b97cf49bSBarry Smith 364433994e6SBarry Smith PetscFunctionBegin; 365b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 366d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 3670f5bd95cSBarry Smith flag = PETSC_FALSE; 368b97cf49bSBarry Smith } 369b97cf49bSBarry Smith 370b97cf49bSBarry Smith /* if the a->i are the same */ 371580bdb30SBarry Smith ierr = PetscArraycmp(a->i,b->i,A->rmap->n+1,&flag);CHKERRQ(ierr); 372b97cf49bSBarry Smith 373b97cf49bSBarry Smith /* if a->j are the same */ 374b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 375b97cf49bSBarry Smith 376b2566f29SBarry Smith ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 3773a40ed3dSBarry Smith PetscFunctionReturn(0); 378b97cf49bSBarry Smith } 379b97cf49bSBarry Smith 38040244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 3816cd91f64SBarry Smith { 382b24ad042SBarry Smith PetscInt i; 3836cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3841a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 3856cd91f64SBarry Smith 3866cd91f64SBarry Smith PetscFunctionBegin; 387d0f46423SBarry Smith *m = A->rmap->n; 3886cd91f64SBarry Smith *ia = a->i; 3896cd91f64SBarry Smith *ja = a->j; 3906cd91f64SBarry Smith *done = PETSC_TRUE; 391d42735a1SBarry Smith if (oshift) { 392d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 393d42735a1SBarry Smith (*ja)[i]++; 394d42735a1SBarry Smith } 395d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 396d42735a1SBarry Smith } 397d42735a1SBarry Smith PetscFunctionReturn(0); 398d42735a1SBarry Smith } 399d42735a1SBarry Smith 40040244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 401d42735a1SBarry Smith { 402b24ad042SBarry Smith PetscInt i; 403d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 4041a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 405d42735a1SBarry Smith 406d42735a1SBarry Smith PetscFunctionBegin; 40770d8d27cSBarry Smith if (ia && a->i != *ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 40870d8d27cSBarry Smith if (ja && a->j != *ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 409d42735a1SBarry Smith if (oshift) { 41070d8d27cSBarry Smith if (!ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument"); 41170d8d27cSBarry Smith if (!ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument"); 412d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 413d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 414d42735a1SBarry Smith (*ja)[i]--; 415d42735a1SBarry Smith } 416d42735a1SBarry Smith } 4176cd91f64SBarry Smith PetscFunctionReturn(0); 4186cd91f64SBarry Smith } 419b97cf49bSBarry Smith 42019fd82e9SBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 42117667f90SBarry Smith { 42217667f90SBarry Smith Mat B; 42317667f90SBarry Smith PetscErrorCode ierr; 42417667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 42517667f90SBarry Smith const PetscInt *rj; 42617667f90SBarry Smith const PetscScalar *ra; 42717667f90SBarry Smith MPI_Comm comm; 42817667f90SBarry Smith 42917667f90SBarry Smith PetscFunctionBegin; 4300298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 4310298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 4320298fd71SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 43317667f90SBarry Smith 43417667f90SBarry Smith /* count the number of nonzeros per row */ 43517667f90SBarry Smith for (i=0; i<m; i++) { 4360298fd71SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 4375ee9ba1cSJed Brown for (j=0; j<len; j++) { 4385ee9ba1cSJed Brown if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 4395ee9ba1cSJed Brown } 44017667f90SBarry Smith nzeros += len; 4410f2063bfSTobin Isaac ierr = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 44217667f90SBarry Smith } 44317667f90SBarry Smith 44417667f90SBarry Smith /* malloc space for nonzeros */ 445854ce69bSBarry Smith ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr); 446854ce69bSBarry Smith ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr); 447854ce69bSBarry Smith ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr); 44817667f90SBarry Smith 44917667f90SBarry Smith nzeros = 0; 45017667f90SBarry Smith ia[0] = 0; 45117667f90SBarry Smith for (i=0; i<m; i++) { 45217667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 45317667f90SBarry Smith cnt = 0; 45417667f90SBarry Smith for (j=0; j<len; j++) { 4555ee9ba1cSJed Brown if (rj[j] != i+rstart) { /* if not diagonal */ 45617667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 45717667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 45817667f90SBarry Smith } 4595ee9ba1cSJed Brown } 46017667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 46117667f90SBarry Smith nzeros += cnt; 46217667f90SBarry Smith ia[i+1] = nzeros; 46317667f90SBarry Smith } 46417667f90SBarry Smith 46517667f90SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 46617667f90SBarry Smith ierr = MatCreate(comm,&B);CHKERRQ(ierr); 46758ec18a5SLisandro Dalcin ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 46817667f90SBarry Smith ierr = MatSetType(B,type);CHKERRQ(ierr); 46917667f90SBarry Smith ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 47017667f90SBarry Smith 471511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 47228be2f97SBarry Smith ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); 47317667f90SBarry Smith } else { 47417667f90SBarry Smith *newmat = B; 47517667f90SBarry Smith } 47617667f90SBarry Smith PetscFunctionReturn(0); 47717667f90SBarry Smith } 47817667f90SBarry Smith 479b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 4800b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 4813eda8832SBarry Smith MatGetRow_MPIAdj, 4823eda8832SBarry Smith MatRestoreRow_MPIAdj, 483b97cf49bSBarry Smith 0, 48497304618SKris Buschelman /* 4*/ 0, 485b97cf49bSBarry Smith 0, 486b97cf49bSBarry Smith 0, 487b97cf49bSBarry Smith 0, 4880b3b1487SBarry Smith 0, 4890b3b1487SBarry Smith 0, 49097304618SKris Buschelman /*10*/ 0, 4910b3b1487SBarry Smith 0, 4920b3b1487SBarry Smith 0, 4930b3b1487SBarry Smith 0, 4940b3b1487SBarry Smith 0, 49597304618SKris Buschelman /*15*/ 0, 4963eda8832SBarry Smith MatEqual_MPIAdj, 4970b3b1487SBarry Smith 0, 4980b3b1487SBarry Smith 0, 4990b3b1487SBarry Smith 0, 50097304618SKris Buschelman /*20*/ 0, 5010b3b1487SBarry Smith 0, 5023eda8832SBarry Smith MatSetOption_MPIAdj, 5030b3b1487SBarry Smith 0, 504d519adbfSMatthew Knepley /*24*/ 0, 5050b3b1487SBarry Smith 0, 5060b3b1487SBarry Smith 0, 5070b3b1487SBarry Smith 0, 5080b3b1487SBarry Smith 0, 509d519adbfSMatthew Knepley /*29*/ 0, 5100b3b1487SBarry Smith 0, 511273d9f13SBarry Smith 0, 512273d9f13SBarry Smith 0, 5130b3b1487SBarry Smith 0, 514d519adbfSMatthew Knepley /*34*/ 0, 5150b3b1487SBarry Smith 0, 5160b3b1487SBarry Smith 0, 5170b3b1487SBarry Smith 0, 5180b3b1487SBarry Smith 0, 519d519adbfSMatthew Knepley /*39*/ 0, 5207dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 5210b3b1487SBarry Smith 0, 5220b3b1487SBarry Smith 0, 5230b3b1487SBarry Smith 0, 524d519adbfSMatthew Knepley /*44*/ 0, 5250b3b1487SBarry Smith 0, 5267d68702bSBarry Smith MatShift_Basic, 5270b3b1487SBarry Smith 0, 5280b3b1487SBarry Smith 0, 529d519adbfSMatthew Knepley /*49*/ 0, 5306cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 531d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 532b97cf49bSBarry Smith 0, 533b97cf49bSBarry Smith 0, 534d519adbfSMatthew Knepley /*54*/ 0, 535b97cf49bSBarry Smith 0, 5360752156aSBarry Smith 0, 5370752156aSBarry Smith 0, 5380b3b1487SBarry Smith 0, 539d519adbfSMatthew Knepley /*59*/ 0, 540b9b97703SBarry Smith MatDestroy_MPIAdj, 541b9b97703SBarry Smith MatView_MPIAdj, 54217667f90SBarry Smith MatConvertFrom_MPIAdj, 54397304618SKris Buschelman 0, 544d519adbfSMatthew Knepley /*64*/ 0, 54597304618SKris Buschelman 0, 54697304618SKris Buschelman 0, 54797304618SKris Buschelman 0, 54897304618SKris Buschelman 0, 549d519adbfSMatthew Knepley /*69*/ 0, 55097304618SKris Buschelman 0, 55197304618SKris Buschelman 0, 55297304618SKris Buschelman 0, 55397304618SKris Buschelman 0, 554d519adbfSMatthew Knepley /*74*/ 0, 55597304618SKris Buschelman 0, 55697304618SKris Buschelman 0, 55797304618SKris Buschelman 0, 55897304618SKris Buschelman 0, 559d519adbfSMatthew Knepley /*79*/ 0, 56097304618SKris Buschelman 0, 56197304618SKris Buschelman 0, 56297304618SKris Buschelman 0, 563865e5f61SKris Buschelman 0, 564d519adbfSMatthew Knepley /*84*/ 0, 565865e5f61SKris Buschelman 0, 566865e5f61SKris Buschelman 0, 567865e5f61SKris Buschelman 0, 568865e5f61SKris Buschelman 0, 569d519adbfSMatthew Knepley /*89*/ 0, 570865e5f61SKris Buschelman 0, 571865e5f61SKris Buschelman 0, 572865e5f61SKris Buschelman 0, 573865e5f61SKris Buschelman 0, 574d519adbfSMatthew Knepley /*94*/ 0, 575865e5f61SKris Buschelman 0, 576865e5f61SKris Buschelman 0, 5773964eb88SJed Brown 0, 5783964eb88SJed Brown 0, 5793964eb88SJed Brown /*99*/ 0, 5803964eb88SJed Brown 0, 5813964eb88SJed Brown 0, 5823964eb88SJed Brown 0, 5833964eb88SJed Brown 0, 5843964eb88SJed Brown /*104*/ 0, 5853964eb88SJed Brown 0, 5863964eb88SJed Brown 0, 5873964eb88SJed Brown 0, 5883964eb88SJed Brown 0, 5893964eb88SJed Brown /*109*/ 0, 5903964eb88SJed Brown 0, 5913964eb88SJed Brown 0, 5923964eb88SJed Brown 0, 5933964eb88SJed Brown 0, 5943964eb88SJed Brown /*114*/ 0, 5953964eb88SJed Brown 0, 5963964eb88SJed Brown 0, 5973964eb88SJed Brown 0, 5983964eb88SJed Brown 0, 5993964eb88SJed Brown /*119*/ 0, 6003964eb88SJed Brown 0, 6013964eb88SJed Brown 0, 6023964eb88SJed Brown 0, 6033964eb88SJed Brown 0, 6043964eb88SJed Brown /*124*/ 0, 6053964eb88SJed Brown 0, 6063964eb88SJed Brown 0, 6073964eb88SJed Brown 0, 6087dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 6093964eb88SJed Brown /*129*/ 0, 6103964eb88SJed Brown 0, 6113964eb88SJed Brown 0, 6123964eb88SJed Brown 0, 6133964eb88SJed Brown 0, 6143964eb88SJed Brown /*134*/ 0, 6153964eb88SJed Brown 0, 6163964eb88SJed Brown 0, 6173964eb88SJed Brown 0, 6183964eb88SJed Brown 0, 6193964eb88SJed Brown /*139*/ 0, 620f9426fe0SMark Adams 0, 6213964eb88SJed Brown 0 6223964eb88SJed Brown }; 623b97cf49bSBarry Smith 624f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 625a23d5eceSKris Buschelman { 626a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 627e895ccc0SFande Kong PetscBool useedgeweights; 628dfbe8321SBarry Smith PetscErrorCode ierr; 6292515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 630b24ad042SBarry Smith PetscInt ii; 631a23d5eceSKris Buschelman #endif 632a23d5eceSKris Buschelman 633a23d5eceSKris Buschelman PetscFunctionBegin; 63426283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 63526283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 636e895ccc0SFande Kong if (values) useedgeweights = PETSC_TRUE; else useedgeweights = PETSC_FALSE; 637e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */ 638e895ccc0SFande Kong ierr = MPIU_Allreduce((int*)&useedgeweights,(int*)&b->useedgeweights,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)B));CHKERRQ(ierr); 63958ec18a5SLisandro Dalcin 6402515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 641e32f2f54SBarry Smith if (i[0] != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]); 642d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 643e02043d6SBarry Smith if (i[ii] < 0 || i[ii] < i[ii-1]) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i[%D]=%D index is out of range: i[%D]=%D",ii,i[ii],ii-1,i[ii-1]); 644a23d5eceSKris Buschelman } 645d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 646e02043d6SBarry Smith if (j[ii] < 0 || j[ii] >= B->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]); 647a23d5eceSKris Buschelman } 648a23d5eceSKris Buschelman #endif 64958ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 650a23d5eceSKris Buschelman 651a23d5eceSKris Buschelman b->j = j; 652a23d5eceSKris Buschelman b->i = i; 653a23d5eceSKris Buschelman b->values = values; 654a23d5eceSKris Buschelman 655d0f46423SBarry Smith b->nz = i[B->rmap->n]; 656a23d5eceSKris Buschelman b->diag = 0; 657a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 658a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 659a23d5eceSKris Buschelman 660a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 661a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 662a23d5eceSKris Buschelman PetscFunctionReturn(0); 663a23d5eceSKris Buschelman } 664b97cf49bSBarry Smith 665f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 6669a3665c6SJed Brown { 6679a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 6689a3665c6SJed Brown PetscErrorCode ierr; 6699a3665c6SJed Brown const PetscInt *ranges; 6709a3665c6SJed Brown MPI_Comm acomm,bcomm; 6719a3665c6SJed Brown MPI_Group agroup,bgroup; 6729a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 6739a3665c6SJed Brown 6749a3665c6SJed Brown PetscFunctionBegin; 6750298fd71SBarry Smith *B = NULL; 676ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr); 6779a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&size);CHKERRQ(ierr); 6789a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr); 6799a3665c6SJed Brown ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 6809a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6819a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 6829a3665c6SJed Brown } 6839a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 6849a3665c6SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 6859a3665c6SJed Brown *B = A; 6869a3665c6SJed Brown PetscFunctionReturn(0); 6879a3665c6SJed Brown } 6889a3665c6SJed Brown 689785e854fSJed Brown ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr); 6909a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6919a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 6929a3665c6SJed Brown } 6939a3665c6SJed Brown ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr); 6949a3665c6SJed Brown ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr); 6959a3665c6SJed Brown ierr = PetscFree(ranks);CHKERRQ(ierr); 6969a3665c6SJed Brown ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr); 6979a3665c6SJed Brown ierr = MPI_Group_free(&agroup);CHKERRQ(ierr); 6989a3665c6SJed Brown ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr); 6999a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 7009a3665c6SJed Brown PetscInt m,N; 7019a3665c6SJed Brown Mat_MPIAdj *b; 7020298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 7030298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 7049a3665c6SJed Brown ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr); 7059a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 7069a3665c6SJed Brown b->freeaij = PETSC_FALSE; 7079a3665c6SJed Brown ierr = MPI_Comm_free(&bcomm);CHKERRQ(ierr); 7089a3665c6SJed Brown } 7099a3665c6SJed Brown PetscFunctionReturn(0); 7109a3665c6SJed Brown } 7119a3665c6SJed Brown 712b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B) 713b44e7856SBarry Smith { 714b44e7856SBarry Smith PetscErrorCode ierr; 715b44e7856SBarry Smith PetscInt M,N,*II,*J,NZ,nz,m,nzstart,i; 716b44e7856SBarry Smith PetscInt *Values = NULL; 717b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 718b44e7856SBarry Smith PetscMPIInt mnz,mm,*allnz,*allm,size,*dispnz,*dispm; 719b44e7856SBarry Smith 720b44e7856SBarry Smith PetscFunctionBegin; 721b44e7856SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 722b44e7856SBarry Smith ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 723b44e7856SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 724b44e7856SBarry Smith nz = adj->nz; 725d4e69552SBarry Smith if (adj->i[m] != nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nz %D not correct i[m] %d",nz,adj->i[m]); 726b44e7856SBarry Smith ierr = MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 727d4e69552SBarry Smith 728b44e7856SBarry Smith ierr = PetscMPIIntCast(nz,&mnz);CHKERRQ(ierr); 729071fcb05SBarry Smith ierr = PetscMalloc2(size,&allnz,size,&dispnz);CHKERRQ(ierr); 730d4e69552SBarry Smith ierr = MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 731b44e7856SBarry Smith dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1]; 732b44e7856SBarry Smith if (adj->values) { 733b44e7856SBarry Smith ierr = PetscMalloc1(NZ,&Values);CHKERRQ(ierr); 734b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 735b44e7856SBarry Smith } 736b44e7856SBarry Smith ierr = PetscMalloc1(NZ,&J);CHKERRQ(ierr); 737b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 738b44e7856SBarry Smith ierr = PetscFree2(allnz,dispnz);CHKERRQ(ierr); 739b44e7856SBarry Smith ierr = MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 740b44e7856SBarry Smith nzstart -= nz; 741b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 742b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] += nzstart; 743b44e7856SBarry Smith ierr = PetscMalloc1(M+1,&II);CHKERRQ(ierr); 744b44e7856SBarry Smith ierr = PetscMPIIntCast(m,&mm);CHKERRQ(ierr); 745b44e7856SBarry Smith ierr = PetscMalloc2(size,&allm,size,&dispm);CHKERRQ(ierr); 746d4e69552SBarry Smith ierr = MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 747b44e7856SBarry Smith dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1]; 748b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 749b44e7856SBarry Smith ierr = PetscFree2(allm,dispm);CHKERRQ(ierr); 750b44e7856SBarry Smith II[M] = NZ; 751b44e7856SBarry Smith /* shift the i[] values back */ 752b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] -= nzstart; 753b44e7856SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B);CHKERRQ(ierr); 754b44e7856SBarry Smith PetscFunctionReturn(0); 755b44e7856SBarry Smith } 756b44e7856SBarry Smith 7579a3665c6SJed Brown /*@ 7589a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 7599a3665c6SJed Brown 7609a3665c6SJed Brown Collective 7619a3665c6SJed Brown 7629a3665c6SJed Brown Input Arguments: 7639a3665c6SJed Brown . A - original MPIAdj matrix 7649a3665c6SJed Brown 7659a3665c6SJed Brown Output Arguments: 7660298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 7679a3665c6SJed Brown 7689a3665c6SJed Brown Level: developer 7699a3665c6SJed Brown 7709a3665c6SJed Brown Note: 7719a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 7729a3665c6SJed Brown 7739a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 7749a3665c6SJed Brown 7759a3665c6SJed Brown .seealso: MatCreateMPIAdj() 7769a3665c6SJed Brown @*/ 7779a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 7789a3665c6SJed Brown { 7799a3665c6SJed Brown PetscErrorCode ierr; 7809a3665c6SJed Brown 7819a3665c6SJed Brown PetscFunctionBegin; 7829a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 7839a3665c6SJed Brown ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 7849a3665c6SJed Brown PetscFunctionReturn(0); 7859a3665c6SJed Brown } 7869a3665c6SJed Brown 7870bad9183SKris Buschelman /*MC 788fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 7890bad9183SKris Buschelman intended for use constructing orderings and partitionings. 7900bad9183SKris Buschelman 7910bad9183SKris Buschelman Level: beginner 7920bad9183SKris Buschelman 7930bad9183SKris Buschelman .seealso: MatCreateMPIAdj 7940bad9183SKris Buschelman M*/ 7950bad9183SKris Buschelman 7968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 797273d9f13SBarry Smith { 798273d9f13SBarry Smith Mat_MPIAdj *b; 7996849ba73SBarry Smith PetscErrorCode ierr; 800273d9f13SBarry Smith 801273d9f13SBarry Smith PetscFunctionBegin; 802b00a9115SJed Brown ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 803b0a32e0cSBarry Smith B->data = (void*)b; 804273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 805273d9f13SBarry Smith B->assembled = PETSC_FALSE; 806273d9f13SBarry Smith 807bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 808bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr); 809b44e7856SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj);CHKERRQ(ierr); 81017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 811273d9f13SBarry Smith PetscFunctionReturn(0); 812273d9f13SBarry Smith } 813273d9f13SBarry Smith 814a23d5eceSKris Buschelman /*@C 815b44e7856SBarry Smith MatMPIAdjToSeq - Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential preconditioners) 816b44e7856SBarry Smith 817d083f849SBarry Smith Logically Collective 818b44e7856SBarry Smith 819b44e7856SBarry Smith Input Parameter: 820b44e7856SBarry Smith . A - the matrix 821b44e7856SBarry Smith 822b44e7856SBarry Smith Output Parameter: 823b44e7856SBarry Smith . B - the same matrix on all processes 824b44e7856SBarry Smith 825b44e7856SBarry Smith Level: intermediate 826b44e7856SBarry Smith 827b44e7856SBarry Smith .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 828b44e7856SBarry Smith @*/ 829b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq(Mat A,Mat *B) 830b44e7856SBarry Smith { 831b44e7856SBarry Smith PetscErrorCode ierr; 832b44e7856SBarry Smith 833b44e7856SBarry Smith PetscFunctionBegin; 834b44e7856SBarry Smith ierr = PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 835b44e7856SBarry Smith PetscFunctionReturn(0); 836b44e7856SBarry Smith } 837b44e7856SBarry Smith 838b44e7856SBarry Smith /*@C 839a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 840a23d5eceSKris Buschelman 841d083f849SBarry Smith Logically Collective 842a23d5eceSKris Buschelman 843a23d5eceSKris Buschelman Input Parameters: 844a23d5eceSKris Buschelman + A - the matrix 845a23d5eceSKris Buschelman . i - the indices into j for the start of each row 846a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 847a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 848a23d5eceSKris Buschelman - values - [optional] edge weights 849a23d5eceSKris Buschelman 850a23d5eceSKris Buschelman Level: intermediate 851a23d5eceSKris Buschelman 852a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 853a23d5eceSKris Buschelman @*/ 8547087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 855273d9f13SBarry Smith { 8564ac538c5SBarry Smith PetscErrorCode ierr; 857273d9f13SBarry Smith 858273d9f13SBarry Smith PetscFunctionBegin; 8594ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 860273d9f13SBarry Smith PetscFunctionReturn(0); 861273d9f13SBarry Smith } 862273d9f13SBarry Smith 863b97cf49bSBarry Smith /*@C 8643eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 865b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 866b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 867b97cf49bSBarry Smith 868d083f849SBarry Smith Collective 869ef5ee4d1SLois Curfman McInnes 870b97cf49bSBarry Smith Input Parameters: 871c2e958feSBarry Smith + comm - MPI communicator 8720752156aSBarry Smith . m - number of local rows 87358ec18a5SLisandro Dalcin . N - number of global columns 874b97cf49bSBarry Smith . i - the indices into j for the start of each row 875329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 876ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 877329f5518SBarry Smith - values -[optional] edge weights 878b97cf49bSBarry Smith 879b97cf49bSBarry Smith Output Parameter: 880b97cf49bSBarry Smith . A - the matrix 881b97cf49bSBarry Smith 88215091d37SBarry Smith Level: intermediate 88315091d37SBarry Smith 88495452b02SPatrick Sanan Notes: 88595452b02SPatrick Sanan This matrix object does not support most matrix operations, include 8864bc6d8c0SBarry Smith MatSetValues(). 887c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 888ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 8891198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 890327e1a73SBarry Smith Should not include the matrix diagonals. 891b97cf49bSBarry Smith 892e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 893ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 894ededeb1aSvictorle 895ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 896b97cf49bSBarry Smith 897e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 898b97cf49bSBarry Smith @*/ 8997087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 900b97cf49bSBarry Smith { 901dfbe8321SBarry Smith PetscErrorCode ierr; 902b97cf49bSBarry Smith 903433994e6SBarry Smith PetscFunctionBegin; 904f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 90558ec18a5SLisandro Dalcin ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 906273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 907273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 9083a40ed3dSBarry Smith PetscFunctionReturn(0); 909b97cf49bSBarry Smith } 910b97cf49bSBarry Smith 911c2e958feSBarry Smith 912433994e6SBarry Smith 913433994e6SBarry Smith 914433994e6SBarry Smith 915433994e6SBarry Smith 916433994e6SBarry Smith 917