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 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); 51ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv,MPI_REPLACE);CHKERRQ(ierr); 52ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv,MPI_REPLACE);CHKERRQ(ierr); 53ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv,MPI_REPLACE);CHKERRQ(ierr); 54ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv,MPI_REPLACE);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); 84ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv,MPI_REPLACE);CHKERRQ(ierr); 85ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv,MPI_REPLACE);CHKERRQ(ierr); 86e895ccc0SFande Kong if (a->useedgeweights) { 87ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv,MPI_REPLACE);CHKERRQ(ierr); 88ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv,MPI_REPLACE);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) { 126f4259b30SLisandro Dalcin if (a->useedgeweights) *sadj_values = svalues; else *sadj_values=NULL; 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); 161ffc4695bSBarry Smith ierr = MPI_Comm_compare(scomm_row,scomm_col,&issame);CHKERRMPI(ierr); 1622c71b3e2SJacob Faibussowitsch PetscCheckFalse(issame != MPI_IDENT,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row index set must have the same comm as the col index set"); 163ffc4695bSBarry Smith ierr = MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame);CHKERRMPI(ierr); 1642c71b3e2SJacob Faibussowitsch PetscCheckFalse(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; 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); 1852c71b3e2SJacob Faibussowitsch PetscCheckFalse(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) { 19040244768SBarry Smith ierr = MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr); 19140244768SBarry Smith } else { 19240244768SBarry Smith Mat sadj = *(submat[i]); 19340244768SBarry Smith Mat_MPIAdj *sa = (Mat_MPIAdj*)((sadj)->data); 19440244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)sadj,&scomm_mat);CHKERRQ(ierr); 195ffc4695bSBarry Smith ierr = MPI_Comm_compare(scomm_row,scomm_mat,&issame);CHKERRMPI(ierr); 1962c71b3e2SJacob Faibussowitsch PetscCheckFalse(issame != MPI_IDENT,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"submatrix must have the same comm as the col index set"); 197580bdb30SBarry Smith ierr = PetscArraycpy(sa->i,sxadj,irow_n+1);CHKERRQ(ierr); 198580bdb30SBarry Smith ierr = PetscArraycpy(sa->j,sadjncy,sxadj[irow_n]);CHKERRQ(ierr); 199580bdb30SBarry Smith if (svalues) {ierr = PetscArraycpy(sa->values,svalues,sxadj[irow_n]);CHKERRQ(ierr);} 20040244768SBarry Smith ierr = PetscFree(sxadj);CHKERRQ(ierr); 20140244768SBarry Smith ierr = PetscFree(sadjncy);CHKERRQ(ierr); 20240244768SBarry Smith if (svalues) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 20340244768SBarry Smith } 20440244768SBarry Smith } 20540244768SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 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 PetscErrorCode ierr; 21240244768SBarry Smith /*get sub-matrices across a sub communicator */ 21340244768SBarry Smith PetscFunctionBegin; 2147dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat);CHKERRQ(ierr); 21540244768SBarry Smith PetscFunctionReturn(0); 21640244768SBarry Smith } 21740244768SBarry Smith 2187dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 21940244768SBarry Smith { 22040244768SBarry Smith PetscErrorCode ierr; 22140244768SBarry Smith 22240244768SBarry Smith PetscFunctionBegin; 22340244768SBarry Smith /*get sub-matrices based on PETSC_COMM_SELF */ 2247dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat);CHKERRQ(ierr); 22540244768SBarry Smith PetscFunctionReturn(0); 22640244768SBarry Smith } 22740244768SBarry Smith 22840244768SBarry Smith static PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 229b97cf49bSBarry Smith { 2303eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 231dfbe8321SBarry Smith PetscErrorCode ierr; 232d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 2332dcb1b2aSMatthew Knepley const char *name; 234ce1411ecSBarry Smith PetscViewerFormat format; 235b97cf49bSBarry Smith 236433994e6SBarry Smith PetscFunctionBegin; 2373a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 238b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 239fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 2403a40ed3dSBarry Smith PetscFunctionReturn(0); 2412c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(format == PETSC_VIEWER_ASCII_MATLAB,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported"); 2426c4ed002SBarry Smith else { 243d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 244c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 245b97cf49bSBarry Smith for (i=0; i<m; i++) { 246c0aa6a63SJacob Faibussowitsch ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %" PetscInt_FMT ":",i+A->rmap->rstart);CHKERRQ(ierr); 247b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 2480170919eSFande Kong if (a->values) { 249c0aa6a63SJacob Faibussowitsch ierr = PetscViewerASCIISynchronizedPrintf(viewer," (%" PetscInt_FMT ", %" PetscInt_FMT ") ",a->j[j], a->values[j]);CHKERRQ(ierr); 2500170919eSFande Kong } else { 251c0aa6a63SJacob Faibussowitsch ierr = PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT " ",a->j[j]);CHKERRQ(ierr); 2520752156aSBarry Smith } 2530170919eSFande Kong } 254b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 255b97cf49bSBarry Smith } 256d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 257b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 258c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2597b23a99aSBarry Smith } 2603a40ed3dSBarry Smith PetscFunctionReturn(0); 261b97cf49bSBarry Smith } 262b97cf49bSBarry Smith 26340244768SBarry Smith static PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 264b97cf49bSBarry Smith { 265dfbe8321SBarry Smith PetscErrorCode ierr; 266ace3abfcSBarry Smith PetscBool iascii; 267b97cf49bSBarry Smith 268433994e6SBarry Smith PetscFunctionBegin; 269251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 27032077d6dSBarry Smith if (iascii) { 2713eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 272b97cf49bSBarry Smith } 2733a40ed3dSBarry Smith PetscFunctionReturn(0); 274b97cf49bSBarry Smith } 275b97cf49bSBarry Smith 27640244768SBarry Smith static PetscErrorCode MatDestroy_MPIAdj(Mat mat) 277b97cf49bSBarry Smith { 2783eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 279dfbe8321SBarry Smith PetscErrorCode ierr; 28094d884c6SBarry Smith 28194d884c6SBarry Smith PetscFunctionBegin; 282aa482453SBarry Smith #if defined(PETSC_USE_LOG) 283c0aa6a63SJacob Faibussowitsch PetscLogObjectState((PetscObject)mat,"Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT,mat->rmap->n,mat->cmap->n,a->nz); 284b97cf49bSBarry Smith #endif 28505b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 2860400557aSBarry Smith if (a->freeaij) { 28714196391SBarry Smith if (a->freeaijwithfree) { 28814196391SBarry Smith if (a->i) free(a->i); 28914196391SBarry Smith if (a->j) free(a->j); 29014196391SBarry Smith } else { 291606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 292606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 29305b42c5fSBarry Smith ierr = PetscFree(a->values);CHKERRQ(ierr); 29414196391SBarry Smith } 2950400557aSBarry Smith } 2962b1d8763SJed Brown ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 297bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 298f4259b30SLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)mat,NULL);CHKERRQ(ierr); 299bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr); 300bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr); 3013a40ed3dSBarry Smith PetscFunctionReturn(0); 302b97cf49bSBarry Smith } 303b97cf49bSBarry Smith 30440244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 305b97cf49bSBarry Smith { 3063eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 30763ba0a88SBarry Smith PetscErrorCode ierr; 308b97cf49bSBarry Smith 309433994e6SBarry Smith PetscFunctionBegin; 31012c028f9SKris Buschelman switch (op) { 31177e54ba9SKris Buschelman case MAT_SYMMETRIC: 31212c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3139a4540c5SBarry Smith case MAT_HERMITIAN: 3144e0d8c25SBarry Smith a->symmetric = flg; 3159a4540c5SBarry Smith break; 3169a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 3179a4540c5SBarry Smith break; 31812c028f9SKris Buschelman default: 3197d3de750SJacob Faibussowitsch ierr = PetscInfo(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 32012c028f9SKris Buschelman break; 321b97cf49bSBarry Smith } 3223a40ed3dSBarry Smith PetscFunctionReturn(0); 323b97cf49bSBarry Smith } 324b97cf49bSBarry Smith 32540244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 326b97cf49bSBarry Smith { 3273eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3282b1d8763SJed Brown PetscErrorCode ierr; 329b97cf49bSBarry Smith 330433994e6SBarry Smith PetscFunctionBegin; 331d0f46423SBarry Smith row -= A->rmap->rstart; 3322c71b3e2SJacob Faibussowitsch PetscCheckFalse(row < 0 || row >= A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 333b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 3342b1d8763SJed Brown if (v) { 3352b1d8763SJed Brown PetscInt j; 3362b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 3372b1d8763SJed Brown ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 3382b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz); 339785e854fSJed Brown ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr); 340b97cf49bSBarry Smith } 34191f8cafeSFande Kong for (j=0; j<*nz; j++) { 34291f8cafeSFande Kong a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0; 34391f8cafeSFande Kong } 3442b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 345b97cf49bSBarry Smith } 34692b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 3473a40ed3dSBarry Smith PetscFunctionReturn(0); 348b97cf49bSBarry Smith } 349b97cf49bSBarry Smith 35040244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 351b97cf49bSBarry Smith { 352433994e6SBarry Smith PetscFunctionBegin; 3533a40ed3dSBarry Smith PetscFunctionReturn(0); 354b97cf49bSBarry Smith } 355b97cf49bSBarry Smith 35640244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 357b97cf49bSBarry Smith { 3583eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 359dfbe8321SBarry Smith PetscErrorCode ierr; 360ace3abfcSBarry Smith PetscBool flag; 361b97cf49bSBarry Smith 362433994e6SBarry Smith PetscFunctionBegin; 363b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 364d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 3650f5bd95cSBarry Smith flag = PETSC_FALSE; 366b97cf49bSBarry Smith } 367b97cf49bSBarry Smith 368b97cf49bSBarry Smith /* if the a->i are the same */ 369580bdb30SBarry Smith ierr = PetscArraycmp(a->i,b->i,A->rmap->n+1,&flag);CHKERRQ(ierr); 370b97cf49bSBarry Smith 371b97cf49bSBarry Smith /* if a->j are the same */ 372b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 373b97cf49bSBarry Smith 374820f2d46SBarry Smith ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 3753a40ed3dSBarry Smith PetscFunctionReturn(0); 376b97cf49bSBarry Smith } 377b97cf49bSBarry Smith 37840244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 3796cd91f64SBarry Smith { 380b24ad042SBarry Smith PetscInt i; 3816cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3821a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 3836cd91f64SBarry Smith 3846cd91f64SBarry Smith PetscFunctionBegin; 385d0f46423SBarry Smith *m = A->rmap->n; 3866cd91f64SBarry Smith *ia = a->i; 3876cd91f64SBarry Smith *ja = a->j; 3886cd91f64SBarry Smith *done = PETSC_TRUE; 389d42735a1SBarry Smith if (oshift) { 390d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 391d42735a1SBarry Smith (*ja)[i]++; 392d42735a1SBarry Smith } 393d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 394d42735a1SBarry Smith } 395d42735a1SBarry Smith PetscFunctionReturn(0); 396d42735a1SBarry Smith } 397d42735a1SBarry Smith 39840244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 399d42735a1SBarry Smith { 400b24ad042SBarry Smith PetscInt i; 401d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 4021a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 403d42735a1SBarry Smith 404d42735a1SBarry Smith PetscFunctionBegin; 4052c71b3e2SJacob Faibussowitsch PetscCheckFalse(ia && a->i != *ia,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 4062c71b3e2SJacob Faibussowitsch PetscCheckFalse(ja && a->j != *ja,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 407d42735a1SBarry Smith if (oshift) { 4082c71b3e2SJacob Faibussowitsch PetscCheckFalse(!ia,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument"); 4092c71b3e2SJacob Faibussowitsch PetscCheckFalse(!ja,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument"); 410d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 411d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 412d42735a1SBarry Smith (*ja)[i]--; 413d42735a1SBarry Smith } 414d42735a1SBarry Smith } 4156cd91f64SBarry Smith PetscFunctionReturn(0); 4166cd91f64SBarry Smith } 417b97cf49bSBarry Smith 41819fd82e9SBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 41917667f90SBarry Smith { 42017667f90SBarry Smith Mat B; 42117667f90SBarry Smith PetscErrorCode ierr; 42217667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 42317667f90SBarry Smith const PetscInt *rj; 42417667f90SBarry Smith const PetscScalar *ra; 42517667f90SBarry Smith MPI_Comm comm; 42617667f90SBarry Smith 42717667f90SBarry Smith PetscFunctionBegin; 4280298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 4290298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 4300298fd71SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 43117667f90SBarry Smith 43217667f90SBarry Smith /* count the number of nonzeros per row */ 43317667f90SBarry Smith for (i=0; i<m; i++) { 4340298fd71SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 4355ee9ba1cSJed Brown for (j=0; j<len; j++) { 4365ee9ba1cSJed Brown if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 4375ee9ba1cSJed Brown } 43817667f90SBarry Smith nzeros += len; 4390f2063bfSTobin Isaac ierr = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 44017667f90SBarry Smith } 44117667f90SBarry Smith 44217667f90SBarry Smith /* malloc space for nonzeros */ 443854ce69bSBarry Smith ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr); 444854ce69bSBarry Smith ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr); 445854ce69bSBarry Smith ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr); 44617667f90SBarry Smith 44717667f90SBarry Smith nzeros = 0; 44817667f90SBarry Smith ia[0] = 0; 44917667f90SBarry Smith for (i=0; i<m; i++) { 45017667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 45117667f90SBarry Smith cnt = 0; 45217667f90SBarry Smith for (j=0; j<len; j++) { 4535ee9ba1cSJed Brown if (rj[j] != i+rstart) { /* if not diagonal */ 45417667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 45517667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 45617667f90SBarry Smith } 4575ee9ba1cSJed Brown } 45817667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 45917667f90SBarry Smith nzeros += cnt; 46017667f90SBarry Smith ia[i+1] = nzeros; 46117667f90SBarry Smith } 46217667f90SBarry Smith 46317667f90SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 46417667f90SBarry Smith ierr = MatCreate(comm,&B);CHKERRQ(ierr); 46558ec18a5SLisandro Dalcin ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 46617667f90SBarry Smith ierr = MatSetType(B,type);CHKERRQ(ierr); 46717667f90SBarry Smith ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 46817667f90SBarry Smith 469511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 47028be2f97SBarry Smith ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); 47117667f90SBarry Smith } else { 47217667f90SBarry Smith *newmat = B; 47317667f90SBarry Smith } 47417667f90SBarry Smith PetscFunctionReturn(0); 47517667f90SBarry Smith } 47617667f90SBarry Smith 477b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 478f4259b30SLisandro Dalcin static struct _MatOps MatOps_Values = {NULL, 4793eda8832SBarry Smith MatGetRow_MPIAdj, 4803eda8832SBarry Smith MatRestoreRow_MPIAdj, 481f4259b30SLisandro Dalcin NULL, 482f4259b30SLisandro Dalcin /* 4*/ NULL, 483f4259b30SLisandro Dalcin NULL, 484f4259b30SLisandro Dalcin NULL, 485f4259b30SLisandro Dalcin NULL, 486f4259b30SLisandro Dalcin NULL, 487f4259b30SLisandro Dalcin NULL, 488f4259b30SLisandro Dalcin /*10*/ NULL, 489f4259b30SLisandro Dalcin NULL, 490f4259b30SLisandro Dalcin NULL, 491f4259b30SLisandro Dalcin NULL, 492f4259b30SLisandro Dalcin NULL, 493f4259b30SLisandro Dalcin /*15*/ NULL, 4943eda8832SBarry Smith MatEqual_MPIAdj, 495f4259b30SLisandro Dalcin NULL, 496f4259b30SLisandro Dalcin NULL, 497f4259b30SLisandro Dalcin NULL, 498f4259b30SLisandro Dalcin /*20*/ NULL, 499f4259b30SLisandro Dalcin NULL, 5003eda8832SBarry Smith MatSetOption_MPIAdj, 501f4259b30SLisandro Dalcin NULL, 502f4259b30SLisandro Dalcin /*24*/ NULL, 503f4259b30SLisandro Dalcin NULL, 504f4259b30SLisandro Dalcin NULL, 505f4259b30SLisandro Dalcin NULL, 506f4259b30SLisandro Dalcin NULL, 507f4259b30SLisandro Dalcin /*29*/ NULL, 508f4259b30SLisandro Dalcin NULL, 509f4259b30SLisandro Dalcin NULL, 510f4259b30SLisandro Dalcin NULL, 511f4259b30SLisandro Dalcin NULL, 512f4259b30SLisandro Dalcin /*34*/ NULL, 513f4259b30SLisandro Dalcin NULL, 514f4259b30SLisandro Dalcin NULL, 515f4259b30SLisandro Dalcin NULL, 516f4259b30SLisandro Dalcin NULL, 517f4259b30SLisandro Dalcin /*39*/ NULL, 5187dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 519f4259b30SLisandro Dalcin NULL, 520f4259b30SLisandro Dalcin NULL, 521f4259b30SLisandro Dalcin NULL, 522f4259b30SLisandro Dalcin /*44*/ NULL, 523f4259b30SLisandro Dalcin NULL, 5247d68702bSBarry Smith MatShift_Basic, 525f4259b30SLisandro Dalcin NULL, 526f4259b30SLisandro Dalcin NULL, 527f4259b30SLisandro Dalcin /*49*/ NULL, 5286cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 529d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 530f4259b30SLisandro Dalcin NULL, 531f4259b30SLisandro Dalcin NULL, 532f4259b30SLisandro Dalcin /*54*/ NULL, 533f4259b30SLisandro Dalcin NULL, 534f4259b30SLisandro Dalcin NULL, 535f4259b30SLisandro Dalcin NULL, 536f4259b30SLisandro Dalcin NULL, 537f4259b30SLisandro Dalcin /*59*/ NULL, 538b9b97703SBarry Smith MatDestroy_MPIAdj, 539b9b97703SBarry Smith MatView_MPIAdj, 54017667f90SBarry Smith MatConvertFrom_MPIAdj, 541f4259b30SLisandro Dalcin NULL, 542f4259b30SLisandro Dalcin /*64*/ NULL, 543f4259b30SLisandro Dalcin NULL, 544f4259b30SLisandro Dalcin NULL, 545f4259b30SLisandro Dalcin NULL, 546f4259b30SLisandro Dalcin NULL, 547f4259b30SLisandro Dalcin /*69*/ NULL, 548f4259b30SLisandro Dalcin NULL, 549f4259b30SLisandro Dalcin NULL, 550f4259b30SLisandro Dalcin NULL, 551f4259b30SLisandro Dalcin NULL, 552f4259b30SLisandro Dalcin /*74*/ NULL, 553f4259b30SLisandro Dalcin NULL, 554f4259b30SLisandro Dalcin NULL, 555f4259b30SLisandro Dalcin NULL, 556f4259b30SLisandro Dalcin NULL, 557f4259b30SLisandro Dalcin /*79*/ NULL, 558f4259b30SLisandro Dalcin NULL, 559f4259b30SLisandro Dalcin NULL, 560f4259b30SLisandro Dalcin NULL, 561f4259b30SLisandro Dalcin NULL, 562f4259b30SLisandro Dalcin /*84*/ NULL, 563f4259b30SLisandro Dalcin NULL, 564f4259b30SLisandro Dalcin NULL, 565f4259b30SLisandro Dalcin NULL, 566f4259b30SLisandro Dalcin NULL, 567f4259b30SLisandro Dalcin /*89*/ NULL, 568f4259b30SLisandro Dalcin NULL, 569f4259b30SLisandro Dalcin NULL, 570f4259b30SLisandro Dalcin NULL, 571f4259b30SLisandro Dalcin NULL, 572f4259b30SLisandro Dalcin /*94*/ NULL, 573f4259b30SLisandro Dalcin NULL, 574f4259b30SLisandro Dalcin NULL, 575f4259b30SLisandro Dalcin NULL, 576f4259b30SLisandro Dalcin NULL, 577f4259b30SLisandro Dalcin /*99*/ NULL, 578f4259b30SLisandro Dalcin NULL, 579f4259b30SLisandro Dalcin NULL, 580f4259b30SLisandro Dalcin NULL, 581f4259b30SLisandro Dalcin NULL, 582f4259b30SLisandro Dalcin /*104*/ NULL, 583f4259b30SLisandro Dalcin NULL, 584f4259b30SLisandro Dalcin NULL, 585f4259b30SLisandro Dalcin NULL, 586f4259b30SLisandro Dalcin NULL, 587f4259b30SLisandro Dalcin /*109*/ NULL, 588f4259b30SLisandro Dalcin NULL, 589f4259b30SLisandro Dalcin NULL, 590f4259b30SLisandro Dalcin NULL, 591f4259b30SLisandro Dalcin NULL, 592f4259b30SLisandro Dalcin /*114*/ NULL, 593f4259b30SLisandro Dalcin NULL, 594f4259b30SLisandro Dalcin NULL, 595f4259b30SLisandro Dalcin NULL, 596f4259b30SLisandro Dalcin NULL, 597f4259b30SLisandro Dalcin /*119*/ NULL, 598f4259b30SLisandro Dalcin NULL, 599f4259b30SLisandro Dalcin NULL, 600f4259b30SLisandro Dalcin NULL, 601f4259b30SLisandro Dalcin NULL, 602f4259b30SLisandro Dalcin /*124*/ NULL, 603f4259b30SLisandro Dalcin NULL, 604f4259b30SLisandro Dalcin NULL, 605f4259b30SLisandro Dalcin NULL, 6067dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 607f4259b30SLisandro Dalcin /*129*/ NULL, 608f4259b30SLisandro Dalcin NULL, 609f4259b30SLisandro Dalcin NULL, 610f4259b30SLisandro Dalcin NULL, 611f4259b30SLisandro Dalcin NULL, 612f4259b30SLisandro Dalcin /*134*/ NULL, 613f4259b30SLisandro Dalcin NULL, 614f4259b30SLisandro Dalcin NULL, 615f4259b30SLisandro Dalcin NULL, 616f4259b30SLisandro Dalcin NULL, 617f4259b30SLisandro Dalcin /*139*/ NULL, 618f4259b30SLisandro Dalcin NULL, 619*d70f29a3SPierre Jolivet NULL, 620*d70f29a3SPierre Jolivet NULL, 621*d70f29a3SPierre Jolivet NULL, 622*d70f29a3SPierre Jolivet /*144*/NULL, 623*d70f29a3SPierre Jolivet NULL, 624*d70f29a3SPierre Jolivet NULL, 625f4259b30SLisandro Dalcin NULL 6263964eb88SJed Brown }; 627b97cf49bSBarry Smith 628f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 629a23d5eceSKris Buschelman { 630a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 631e895ccc0SFande Kong PetscBool useedgeweights; 632dfbe8321SBarry Smith PetscErrorCode ierr; 633a23d5eceSKris Buschelman 634a23d5eceSKris Buschelman PetscFunctionBegin; 63526283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 63626283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 637e895ccc0SFande Kong if (values) useedgeweights = PETSC_TRUE; else useedgeweights = PETSC_FALSE; 638e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */ 639820f2d46SBarry Smith ierr = MPIU_Allreduce((int*)&useedgeweights,(int*)&b->useedgeweights,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)B));CHKERRMPI(ierr); 64058ec18a5SLisandro Dalcin 64176bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 64276bd3646SJed Brown PetscInt ii; 64376bd3646SJed Brown 6442c71b3e2SJacob Faibussowitsch PetscCheckFalse(i[0] != 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %" PetscInt_FMT,i[0]); 645d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 6462c71b3e2SJacob Faibussowitsch PetscCheckFalse(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]); 647a23d5eceSKris Buschelman } 648d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 6492c71b3e2SJacob Faibussowitsch PetscCheckFalse(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]); 650a23d5eceSKris Buschelman } 65176bd3646SJed Brown } 65258ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 653a23d5eceSKris Buschelman 654a23d5eceSKris Buschelman b->j = j; 655a23d5eceSKris Buschelman b->i = i; 656a23d5eceSKris Buschelman b->values = values; 657a23d5eceSKris Buschelman 658d0f46423SBarry Smith b->nz = i[B->rmap->n]; 659f4259b30SLisandro Dalcin b->diag = NULL; 660a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 661a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 662a23d5eceSKris Buschelman 663a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 664a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 665a23d5eceSKris Buschelman PetscFunctionReturn(0); 666a23d5eceSKris Buschelman } 667b97cf49bSBarry Smith 668f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 6699a3665c6SJed Brown { 6709a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 6719a3665c6SJed Brown PetscErrorCode ierr; 6729a3665c6SJed Brown const PetscInt *ranges; 6739a3665c6SJed Brown MPI_Comm acomm,bcomm; 6749a3665c6SJed Brown MPI_Group agroup,bgroup; 6759a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 6769a3665c6SJed Brown 6779a3665c6SJed Brown PetscFunctionBegin; 6780298fd71SBarry Smith *B = NULL; 679ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr); 68055b25c41SPierre Jolivet ierr = MPI_Comm_size(acomm,&size);CHKERRMPI(ierr); 68155b25c41SPierre Jolivet ierr = MPI_Comm_size(acomm,&rank);CHKERRMPI(ierr); 6829a3665c6SJed Brown ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 6839a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6849a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 6859a3665c6SJed Brown } 6869a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 6879a3665c6SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 6889a3665c6SJed Brown *B = A; 6899a3665c6SJed Brown PetscFunctionReturn(0); 6909a3665c6SJed Brown } 6919a3665c6SJed Brown 692785e854fSJed Brown ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr); 6939a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6949a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 6959a3665c6SJed Brown } 696ffc4695bSBarry Smith ierr = MPI_Comm_group(acomm,&agroup);CHKERRMPI(ierr); 697ffc4695bSBarry Smith ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRMPI(ierr); 6989a3665c6SJed Brown ierr = PetscFree(ranks);CHKERRQ(ierr); 699ffc4695bSBarry Smith ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRMPI(ierr); 700ffc4695bSBarry Smith ierr = MPI_Group_free(&agroup);CHKERRMPI(ierr); 701ffc4695bSBarry Smith ierr = MPI_Group_free(&bgroup);CHKERRMPI(ierr); 7029a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 7039a3665c6SJed Brown PetscInt m,N; 7049a3665c6SJed Brown Mat_MPIAdj *b; 7050298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 7060298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 7079a3665c6SJed Brown ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr); 7089a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 7099a3665c6SJed Brown b->freeaij = PETSC_FALSE; 71055b25c41SPierre Jolivet ierr = MPI_Comm_free(&bcomm);CHKERRMPI(ierr); 7119a3665c6SJed Brown } 7129a3665c6SJed Brown PetscFunctionReturn(0); 7139a3665c6SJed Brown } 7149a3665c6SJed Brown 715b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B) 716b44e7856SBarry Smith { 717b44e7856SBarry Smith PetscErrorCode ierr; 718b44e7856SBarry Smith PetscInt M,N,*II,*J,NZ,nz,m,nzstart,i; 719b44e7856SBarry Smith PetscInt *Values = NULL; 720b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 721b44e7856SBarry Smith PetscMPIInt mnz,mm,*allnz,*allm,size,*dispnz,*dispm; 722b44e7856SBarry Smith 723b44e7856SBarry Smith PetscFunctionBegin; 724ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 725b44e7856SBarry Smith ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 726b44e7856SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 727b44e7856SBarry Smith nz = adj->nz; 7282c71b3e2SJacob Faibussowitsch PetscCheckFalse(adj->i[m] != nz,PETSC_COMM_SELF,PETSC_ERR_PLIB,"nz %" PetscInt_FMT " not correct i[m] %" PetscInt_FMT,nz,adj->i[m]); 729ffc4695bSBarry Smith ierr = MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 730d4e69552SBarry Smith 731b44e7856SBarry Smith ierr = PetscMPIIntCast(nz,&mnz);CHKERRQ(ierr); 732071fcb05SBarry Smith ierr = PetscMalloc2(size,&allnz,size,&dispnz);CHKERRQ(ierr); 733ffc4695bSBarry Smith ierr = MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 734b44e7856SBarry Smith dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1]; 735b44e7856SBarry Smith if (adj->values) { 736b44e7856SBarry Smith ierr = PetscMalloc1(NZ,&Values);CHKERRQ(ierr); 737ffc4695bSBarry Smith ierr = MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 738b44e7856SBarry Smith } 739b44e7856SBarry Smith ierr = PetscMalloc1(NZ,&J);CHKERRQ(ierr); 740ffc4695bSBarry Smith ierr = MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 741b44e7856SBarry Smith ierr = PetscFree2(allnz,dispnz);CHKERRQ(ierr); 742ffc4695bSBarry Smith ierr = MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 743b44e7856SBarry Smith nzstart -= nz; 744b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 745b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] += nzstart; 746b44e7856SBarry Smith ierr = PetscMalloc1(M+1,&II);CHKERRQ(ierr); 747b44e7856SBarry Smith ierr = PetscMPIIntCast(m,&mm);CHKERRQ(ierr); 748b44e7856SBarry Smith ierr = PetscMalloc2(size,&allm,size,&dispm);CHKERRQ(ierr); 749ffc4695bSBarry Smith ierr = MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 750b44e7856SBarry Smith dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1]; 751ffc4695bSBarry Smith ierr = MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRMPI(ierr); 752b44e7856SBarry Smith ierr = PetscFree2(allm,dispm);CHKERRQ(ierr); 753b44e7856SBarry Smith II[M] = NZ; 754b44e7856SBarry Smith /* shift the i[] values back */ 755b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] -= nzstart; 756b44e7856SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B);CHKERRQ(ierr); 757b44e7856SBarry Smith PetscFunctionReturn(0); 758b44e7856SBarry Smith } 759b44e7856SBarry Smith 7609a3665c6SJed Brown /*@ 7619a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 7629a3665c6SJed Brown 7639a3665c6SJed Brown Collective 7649a3665c6SJed Brown 7654165533cSJose E. Roman Input Parameter: 7669a3665c6SJed Brown . A - original MPIAdj matrix 7679a3665c6SJed Brown 7684165533cSJose E. Roman Output Parameter: 7690298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 7709a3665c6SJed Brown 7719a3665c6SJed Brown Level: developer 7729a3665c6SJed Brown 7739a3665c6SJed Brown Note: 7749a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 7759a3665c6SJed Brown 7769a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 7779a3665c6SJed Brown 7789a3665c6SJed Brown .seealso: MatCreateMPIAdj() 7799a3665c6SJed Brown @*/ 7809a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 7819a3665c6SJed Brown { 7829a3665c6SJed Brown PetscErrorCode ierr; 7839a3665c6SJed Brown 7849a3665c6SJed Brown PetscFunctionBegin; 7859a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 7869a3665c6SJed Brown ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 7879a3665c6SJed Brown PetscFunctionReturn(0); 7889a3665c6SJed Brown } 7899a3665c6SJed Brown 7900bad9183SKris Buschelman /*MC 791fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 7920bad9183SKris Buschelman intended for use constructing orderings and partitionings. 7930bad9183SKris Buschelman 7940bad9183SKris Buschelman Level: beginner 7950bad9183SKris Buschelman 7960bad9183SKris Buschelman .seealso: MatCreateMPIAdj 7970bad9183SKris Buschelman M*/ 7980bad9183SKris Buschelman 7998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 800273d9f13SBarry Smith { 801273d9f13SBarry Smith Mat_MPIAdj *b; 8026849ba73SBarry Smith PetscErrorCode ierr; 803273d9f13SBarry Smith 804273d9f13SBarry Smith PetscFunctionBegin; 805b00a9115SJed Brown ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 806b0a32e0cSBarry Smith B->data = (void*)b; 807273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 808273d9f13SBarry Smith B->assembled = PETSC_FALSE; 809273d9f13SBarry Smith 810bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 811bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr); 812b44e7856SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj);CHKERRQ(ierr); 81317667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 814273d9f13SBarry Smith PetscFunctionReturn(0); 815273d9f13SBarry Smith } 816273d9f13SBarry Smith 817a23d5eceSKris Buschelman /*@C 818b44e7856SBarry Smith MatMPIAdjToSeq - Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential preconditioners) 819b44e7856SBarry Smith 820d083f849SBarry Smith Logically Collective 821b44e7856SBarry Smith 822b44e7856SBarry Smith Input Parameter: 823b44e7856SBarry Smith . A - the matrix 824b44e7856SBarry Smith 825b44e7856SBarry Smith Output Parameter: 826b44e7856SBarry Smith . B - the same matrix on all processes 827b44e7856SBarry Smith 828b44e7856SBarry Smith Level: intermediate 829b44e7856SBarry Smith 830b44e7856SBarry Smith .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 831b44e7856SBarry Smith @*/ 832b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq(Mat A,Mat *B) 833b44e7856SBarry Smith { 834b44e7856SBarry Smith PetscErrorCode ierr; 835b44e7856SBarry Smith 836b44e7856SBarry Smith PetscFunctionBegin; 837b44e7856SBarry Smith ierr = PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 838b44e7856SBarry Smith PetscFunctionReturn(0); 839b44e7856SBarry Smith } 840b44e7856SBarry Smith 841b44e7856SBarry Smith /*@C 842a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 843a23d5eceSKris Buschelman 844d083f849SBarry Smith Logically Collective 845a23d5eceSKris Buschelman 846a23d5eceSKris Buschelman Input Parameters: 847a23d5eceSKris Buschelman + A - the matrix 848a23d5eceSKris Buschelman . i - the indices into j for the start of each row 849a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 850a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 851a23d5eceSKris Buschelman - values - [optional] edge weights 852a23d5eceSKris Buschelman 853a23d5eceSKris Buschelman Level: intermediate 854a23d5eceSKris Buschelman 855a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 856a23d5eceSKris Buschelman @*/ 8577087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 858273d9f13SBarry Smith { 8594ac538c5SBarry Smith PetscErrorCode ierr; 860273d9f13SBarry Smith 861273d9f13SBarry Smith PetscFunctionBegin; 8624ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 863273d9f13SBarry Smith PetscFunctionReturn(0); 864273d9f13SBarry Smith } 865273d9f13SBarry Smith 866b97cf49bSBarry Smith /*@C 8673eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 868b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 869b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 870b97cf49bSBarry Smith 871d083f849SBarry Smith Collective 872ef5ee4d1SLois Curfman McInnes 873b97cf49bSBarry Smith Input Parameters: 874c2e958feSBarry Smith + comm - MPI communicator 8750752156aSBarry Smith . m - number of local rows 87658ec18a5SLisandro Dalcin . N - number of global columns 877b97cf49bSBarry Smith . i - the indices into j for the start of each row 878329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 879ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 880329f5518SBarry Smith - values -[optional] edge weights 881b97cf49bSBarry Smith 882b97cf49bSBarry Smith Output Parameter: 883b97cf49bSBarry Smith . A - the matrix 884b97cf49bSBarry Smith 88515091d37SBarry Smith Level: intermediate 88615091d37SBarry Smith 88795452b02SPatrick Sanan Notes: 88895452b02SPatrick Sanan This matrix object does not support most matrix operations, include 8894bc6d8c0SBarry Smith MatSetValues(). 890c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 891ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 8921198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 893327e1a73SBarry Smith Should not include the matrix diagonals. 894b97cf49bSBarry Smith 895e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 896ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 897ededeb1aSvictorle 898ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 899b97cf49bSBarry Smith 900e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 901b97cf49bSBarry Smith @*/ 9027087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 903b97cf49bSBarry Smith { 904dfbe8321SBarry Smith PetscErrorCode ierr; 905b97cf49bSBarry Smith 906433994e6SBarry Smith PetscFunctionBegin; 907f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 90858ec18a5SLisandro Dalcin ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 909273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 910273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 9113a40ed3dSBarry Smith PetscFunctionReturn(0); 912b97cf49bSBarry Smith } 913