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 { 1340244768SBarry Smith PetscInt nlrows_is,icols_n,i,j,nroots,nleaves,owner,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; 1740244768SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)adj->data; 1840244768SBarry Smith PetscLayout rmap; 1940244768SBarry Smith MPI_Comm comm; 2040244768SBarry Smith PetscSF sf; 2140244768SBarry Smith PetscSFNode *iremote; 2240244768SBarry Smith PetscBool done; 2340244768SBarry Smith PetscErrorCode ierr; 2440244768SBarry Smith 2540244768SBarry Smith PetscFunctionBegin; 2640244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)adj,&comm);CHKERRQ(ierr); 27390e1bf2SBarry Smith ierr = MatGetLayouts(adj,&rmap,NULL);CHKERRQ(ierr); 2840244768SBarry Smith ierr = ISGetLocalSize(irows,&nlrows_is);CHKERRQ(ierr); 2940244768SBarry Smith ierr = ISGetIndices(irows,&irows_indices);CHKERRQ(ierr); 30*071fcb05SBarry Smith ierr = PetscMalloc1(nlrows_is,&iremote);CHKERRQ(ierr); 3140244768SBarry Smith /* construct sf graph*/ 3240244768SBarry Smith nleaves = nlrows_is; 3340244768SBarry Smith for (i=0; i<nlrows_is; i++){ 3440244768SBarry Smith owner = -1; 3540244768SBarry Smith rlocalindex = -1; 3640244768SBarry Smith ierr = PetscLayoutFindOwnerIndex(rmap,irows_indices[i],&owner,&rlocalindex);CHKERRQ(ierr); 3740244768SBarry Smith iremote[i].rank = owner; 3840244768SBarry Smith iremote[i].index = rlocalindex; 3940244768SBarry Smith } 4040244768SBarry Smith ierr = MatGetRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr); 4140244768SBarry Smith ierr = PetscCalloc4(nlrows_mat,&ncols_send,nlrows_is,&xadj_recv,nlrows_is+1,&ncols_recv_offsets,nlrows_is,&ncols_recv);CHKERRQ(ierr); 4240244768SBarry Smith nroots = nlrows_mat; 4340244768SBarry Smith for (i=0; i<nlrows_mat; i++){ 4440244768SBarry Smith ncols_send[i] = xadj[i+1]-xadj[i]; 4540244768SBarry Smith } 4640244768SBarry Smith ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 47390e1bf2SBarry Smith ierr = PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 4840244768SBarry Smith ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 4940244768SBarry Smith ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 5040244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr); 5140244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr); 5240244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr); 5340244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr); 5440244768SBarry Smith ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 5540244768SBarry Smith Ncols_recv =0; 5640244768SBarry Smith for (i=0; i<nlrows_is; i++){ 5740244768SBarry Smith Ncols_recv += ncols_recv[i]; 5840244768SBarry Smith ncols_recv_offsets[i+1] = ncols_recv[i]+ncols_recv_offsets[i]; 5940244768SBarry Smith } 6040244768SBarry Smith Ncols_send = 0; 6140244768SBarry Smith for(i=0; i<nlrows_mat; i++){ 6240244768SBarry Smith Ncols_send += ncols_send[i]; 6340244768SBarry Smith } 6440244768SBarry Smith ierr = PetscCalloc1(Ncols_recv,&iremote);CHKERRQ(ierr); 6540244768SBarry Smith ierr = PetscCalloc1(Ncols_recv,&adjncy_recv);CHKERRQ(ierr); 6640244768SBarry Smith nleaves = Ncols_recv; 6740244768SBarry Smith Ncols_recv = 0; 6840244768SBarry Smith for (i=0; i<nlrows_is; i++){ 6940244768SBarry Smith ierr = PetscLayoutFindOwner(rmap,irows_indices[i],&owner);CHKERRQ(ierr); 7040244768SBarry Smith for (j=0; j<ncols_recv[i]; j++){ 7140244768SBarry Smith iremote[Ncols_recv].rank = owner; 7240244768SBarry Smith iremote[Ncols_recv++].index = xadj_recv[i]+j; 7340244768SBarry Smith } 7440244768SBarry Smith } 7540244768SBarry Smith ierr = ISRestoreIndices(irows,&irows_indices);CHKERRQ(ierr); 7640244768SBarry Smith /*if we need to deal with edge weights ???*/ 77e895ccc0SFande Kong if (a->useedgeweights){ierr = PetscCalloc1(Ncols_recv,&values_recv);CHKERRQ(ierr);} 7840244768SBarry Smith nroots = Ncols_send; 7940244768SBarry Smith ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 80390e1bf2SBarry Smith ierr = PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 8140244768SBarry Smith ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 8240244768SBarry Smith ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 8340244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 8440244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 85e895ccc0SFande Kong if (a->useedgeweights){ 8640244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 8740244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 8840244768SBarry Smith } 8940244768SBarry Smith ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 9040244768SBarry Smith ierr = MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr); 9140244768SBarry Smith ierr = ISGetLocalSize(icols,&icols_n);CHKERRQ(ierr); 9240244768SBarry Smith ierr = ISGetIndices(icols,&icols_indices);CHKERRQ(ierr); 9340244768SBarry Smith rnclos = 0; 9440244768SBarry Smith for (i=0; i<nlrows_is; i++){ 9540244768SBarry Smith for (j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 9640244768SBarry Smith ierr = PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc);CHKERRQ(ierr); 9740244768SBarry Smith if (loc<0){ 9840244768SBarry Smith adjncy_recv[j] = -1; 99e895ccc0SFande Kong if (a->useedgeweights) values_recv[j] = -1; 10040244768SBarry Smith ncols_recv[i]--; 10140244768SBarry Smith } else { 10240244768SBarry Smith rnclos++; 10340244768SBarry Smith } 10440244768SBarry Smith } 10540244768SBarry Smith } 10640244768SBarry Smith ierr = ISRestoreIndices(icols,&icols_indices);CHKERRQ(ierr); 10740244768SBarry Smith ierr = PetscCalloc1(rnclos,&sadjncy);CHKERRQ(ierr); 108e895ccc0SFande Kong if (a->useedgeweights) {ierr = PetscCalloc1(rnclos,&svalues);CHKERRQ(ierr);} 10940244768SBarry Smith ierr = PetscCalloc1(nlrows_is+1,&sxadj);CHKERRQ(ierr); 11040244768SBarry Smith rnclos = 0; 11140244768SBarry Smith for(i=0; i<nlrows_is; i++){ 11240244768SBarry Smith for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 11340244768SBarry Smith if (adjncy_recv[j]<0) continue; 11440244768SBarry Smith sadjncy[rnclos] = adjncy_recv[j]; 115e895ccc0SFande Kong if (a->useedgeweights) svalues[rnclos] = values_recv[j]; 11640244768SBarry Smith rnclos++; 11740244768SBarry Smith } 11840244768SBarry Smith } 11940244768SBarry Smith for(i=0; i<nlrows_is; i++){ 12040244768SBarry Smith sxadj[i+1] = sxadj[i]+ncols_recv[i]; 12140244768SBarry Smith } 12240244768SBarry Smith if (sadj_xadj) { *sadj_xadj = sxadj;} else { ierr = PetscFree(sxadj);CHKERRQ(ierr);} 12340244768SBarry Smith if (sadj_adjncy){ *sadj_adjncy = sadjncy;} else { ierr = PetscFree(sadjncy);CHKERRQ(ierr);} 12440244768SBarry Smith if (sadj_values){ 125e895ccc0SFande Kong if (a->useedgeweights) *sadj_values = svalues; else *sadj_values=0; 12640244768SBarry Smith } else { 127e895ccc0SFande Kong if (a->useedgeweights) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 12840244768SBarry Smith } 12940244768SBarry Smith ierr = PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv);CHKERRQ(ierr); 13040244768SBarry Smith ierr = PetscFree(adjncy_recv);CHKERRQ(ierr); 131e895ccc0SFande Kong if (a->useedgeweights) {ierr = PetscFree(values_recv);CHKERRQ(ierr);} 13240244768SBarry Smith PetscFunctionReturn(0); 13340244768SBarry Smith } 13440244768SBarry Smith 1357dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj_Private(Mat mat,PetscInt n,const IS irow[],const IS icol[],PetscBool subcomm,MatReuse scall,Mat *submat[]) 13640244768SBarry Smith { 13740244768SBarry Smith PetscInt i,irow_n,icol_n,*sxadj,*sadjncy,*svalues; 13840244768SBarry Smith PetscInt *indices,nindx,j,k,loc; 13940244768SBarry Smith PetscMPIInt issame; 14040244768SBarry Smith const PetscInt *irow_indices,*icol_indices; 14140244768SBarry Smith MPI_Comm scomm_row,scomm_col,scomm_mat; 14240244768SBarry Smith PetscErrorCode ierr; 14340244768SBarry Smith 14440244768SBarry Smith PetscFunctionBegin; 14540244768SBarry Smith nindx = 0; 14640244768SBarry Smith /* 14740244768SBarry Smith * Estimate a maximum number for allocating memory 14840244768SBarry Smith */ 14940244768SBarry Smith for(i=0; i<n; i++){ 15040244768SBarry Smith ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr); 15140244768SBarry Smith ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr); 15240244768SBarry Smith nindx = nindx>(irow_n+icol_n)? nindx:(irow_n+icol_n); 15340244768SBarry Smith } 154*071fcb05SBarry Smith ierr = PetscMalloc1(nindx,&indices);CHKERRQ(ierr); 15540244768SBarry Smith /* construct a submat */ 15640244768SBarry Smith for (i=0; i<n; i++){ 15740244768SBarry Smith if (subcomm){ 15840244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)irow[i],&scomm_row);CHKERRQ(ierr); 15940244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)icol[i],&scomm_col);CHKERRQ(ierr); 16040244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,scomm_col,&issame);CHKERRQ(ierr); 16140244768SBarry 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"); 16240244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame);CHKERRQ(ierr); 16340244768SBarry 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"); 16440244768SBarry Smith } else { 16540244768SBarry Smith scomm_row = PETSC_COMM_SELF; 16640244768SBarry Smith } 16740244768SBarry Smith /*get sub-matrix data*/ 16840244768SBarry Smith sxadj=0; sadjncy=0; svalues=0; 1697dae84e0SHong Zhang ierr = MatCreateSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues);CHKERRQ(ierr); 17040244768SBarry Smith ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr); 17140244768SBarry Smith ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr); 17240244768SBarry Smith ierr = ISGetIndices(irow[i],&irow_indices);CHKERRQ(ierr); 173580bdb30SBarry Smith ierr = PetscArraycpy(indices,irow_indices,irow_n);CHKERRQ(ierr); 17440244768SBarry Smith ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr); 17540244768SBarry Smith ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr); 176580bdb30SBarry Smith ierr = PetscArraycpy(indices+irow_n,icol_indices,icol_n);CHKERRQ(ierr); 17740244768SBarry Smith ierr = ISRestoreIndices(icol[i],&icol_indices);CHKERRQ(ierr); 17840244768SBarry Smith nindx = irow_n+icol_n; 17940244768SBarry Smith ierr = PetscSortRemoveDupsInt(&nindx,indices);CHKERRQ(ierr); 18040244768SBarry Smith /* renumber columns */ 18140244768SBarry Smith for (j=0; j<irow_n; j++){ 18240244768SBarry Smith for (k=sxadj[j]; k<sxadj[j+1]; k++){ 18340244768SBarry Smith ierr = PetscFindInt(sadjncy[k],nindx,indices,&loc);CHKERRQ(ierr); 18462795ce3SStefano Zampini #if defined(PETSC_USE_DEBUG) 18562795ce3SStefano Zampini if (loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %D",sadjncy[k]); 18640244768SBarry Smith #endif 18740244768SBarry Smith sadjncy[k] = loc; 18840244768SBarry Smith } 18940244768SBarry Smith } 19040244768SBarry Smith if (scall==MAT_INITIAL_MATRIX){ 19140244768SBarry Smith ierr = MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr); 19240244768SBarry Smith } else { 19340244768SBarry Smith Mat sadj = *(submat[i]); 19440244768SBarry Smith Mat_MPIAdj *sa = (Mat_MPIAdj*)((sadj)->data); 19540244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)sadj,&scomm_mat);CHKERRQ(ierr); 19640244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,scomm_mat,&issame);CHKERRQ(ierr); 19740244768SBarry 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"); 198580bdb30SBarry Smith ierr = PetscArraycpy(sa->i,sxadj,irow_n+1);CHKERRQ(ierr); 199580bdb30SBarry Smith ierr = PetscArraycpy(sa->j,sadjncy,sxadj[irow_n]);CHKERRQ(ierr); 200580bdb30SBarry Smith if (svalues){ierr = PetscArraycpy(sa->values,svalues,sxadj[irow_n]);CHKERRQ(ierr);} 20140244768SBarry Smith ierr = PetscFree(sxadj);CHKERRQ(ierr); 20240244768SBarry Smith ierr = PetscFree(sadjncy);CHKERRQ(ierr); 20340244768SBarry Smith if (svalues) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 20440244768SBarry Smith } 20540244768SBarry Smith } 20640244768SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 20740244768SBarry Smith PetscFunctionReturn(0); 20840244768SBarry Smith } 20940244768SBarry Smith 2107dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatricesMPI_MPIAdj(Mat mat,PetscInt n, const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 21140244768SBarry Smith { 21240244768SBarry Smith PetscErrorCode ierr; 21340244768SBarry Smith /*get sub-matrices across a sub communicator */ 21440244768SBarry Smith PetscFunctionBegin; 2157dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat);CHKERRQ(ierr); 21640244768SBarry Smith PetscFunctionReturn(0); 21740244768SBarry Smith } 21840244768SBarry Smith 2197dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 22040244768SBarry Smith { 22140244768SBarry Smith PetscErrorCode ierr; 22240244768SBarry Smith 22340244768SBarry Smith PetscFunctionBegin; 22440244768SBarry Smith /*get sub-matrices based on PETSC_COMM_SELF */ 2257dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat);CHKERRQ(ierr); 22640244768SBarry Smith PetscFunctionReturn(0); 22740244768SBarry Smith } 22840244768SBarry Smith 22940244768SBarry Smith static PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 230b97cf49bSBarry Smith { 2313eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 232dfbe8321SBarry Smith PetscErrorCode ierr; 233d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 2342dcb1b2aSMatthew Knepley const char *name; 235ce1411ecSBarry Smith PetscViewerFormat format; 236b97cf49bSBarry Smith 237433994e6SBarry Smith PetscFunctionBegin; 2383a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 239b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 240fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 2413a40ed3dSBarry Smith PetscFunctionReturn(0); 2426c4ed002SBarry Smith } else if (format == PETSC_VIEWER_ASCII_MATLAB) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported"); 2436c4ed002SBarry Smith else { 244d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 245c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 246b97cf49bSBarry Smith for (i=0; i<m; i++) { 247d0f46423SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr); 248b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 2490170919eSFande Kong if (a->values) { 2500170919eSFande Kong ierr = PetscViewerASCIISynchronizedPrintf(viewer," (%D, %D) ",a->j[j], a->values[j]);CHKERRQ(ierr); 2510170919eSFande Kong } else { 25277431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 2530752156aSBarry Smith } 2540170919eSFande Kong } 255b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 256b97cf49bSBarry Smith } 257d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 258b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 259c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2607b23a99aSBarry Smith } 2613a40ed3dSBarry Smith PetscFunctionReturn(0); 262b97cf49bSBarry Smith } 263b97cf49bSBarry Smith 26440244768SBarry Smith static PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 265b97cf49bSBarry Smith { 266dfbe8321SBarry Smith PetscErrorCode ierr; 267ace3abfcSBarry Smith PetscBool iascii; 268b97cf49bSBarry Smith 269433994e6SBarry Smith PetscFunctionBegin; 270251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 27132077d6dSBarry Smith if (iascii) { 2723eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 273b97cf49bSBarry Smith } 2743a40ed3dSBarry Smith PetscFunctionReturn(0); 275b97cf49bSBarry Smith } 276b97cf49bSBarry Smith 27740244768SBarry Smith static PetscErrorCode MatDestroy_MPIAdj(Mat mat) 278b97cf49bSBarry Smith { 2793eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 280dfbe8321SBarry Smith PetscErrorCode ierr; 28194d884c6SBarry Smith 28294d884c6SBarry Smith PetscFunctionBegin; 283aa482453SBarry Smith #if defined(PETSC_USE_LOG) 284d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz); 285b97cf49bSBarry Smith #endif 28605b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 2870400557aSBarry Smith if (a->freeaij) { 28814196391SBarry Smith if (a->freeaijwithfree) { 28914196391SBarry Smith if (a->i) free(a->i); 29014196391SBarry Smith if (a->j) free(a->j); 29114196391SBarry Smith } else { 292606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 293606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 29405b42c5fSBarry Smith ierr = PetscFree(a->values);CHKERRQ(ierr); 29514196391SBarry Smith } 2960400557aSBarry Smith } 2972b1d8763SJed Brown ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 298bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 299dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 300bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr); 301bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr); 3023a40ed3dSBarry Smith PetscFunctionReturn(0); 303b97cf49bSBarry Smith } 304b97cf49bSBarry Smith 30540244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 306b97cf49bSBarry Smith { 3073eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 30863ba0a88SBarry Smith PetscErrorCode ierr; 309b97cf49bSBarry Smith 310433994e6SBarry Smith PetscFunctionBegin; 31112c028f9SKris Buschelman switch (op) { 31277e54ba9SKris Buschelman case MAT_SYMMETRIC: 31312c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3149a4540c5SBarry Smith case MAT_HERMITIAN: 3154e0d8c25SBarry Smith a->symmetric = flg; 3169a4540c5SBarry Smith break; 3179a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 3189a4540c5SBarry Smith break; 31912c028f9SKris Buschelman default: 320290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 32112c028f9SKris Buschelman break; 322b97cf49bSBarry Smith } 3233a40ed3dSBarry Smith PetscFunctionReturn(0); 324b97cf49bSBarry Smith } 325b97cf49bSBarry Smith 32640244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 327b97cf49bSBarry Smith { 3283eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3292b1d8763SJed Brown PetscErrorCode ierr; 330b97cf49bSBarry Smith 331433994e6SBarry Smith PetscFunctionBegin; 332d0f46423SBarry Smith row -= A->rmap->rstart; 333e32f2f54SBarry Smith if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 334b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 3352b1d8763SJed Brown if (v) { 3362b1d8763SJed Brown PetscInt j; 3372b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 3382b1d8763SJed Brown ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 3392b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz); 340785e854fSJed Brown ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr); 341b97cf49bSBarry Smith } 34291f8cafeSFande Kong for (j=0; j<*nz; j++){ 34391f8cafeSFande Kong a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0; 34491f8cafeSFande Kong } 3452b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 346b97cf49bSBarry Smith } 34792b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 3483a40ed3dSBarry Smith PetscFunctionReturn(0); 349b97cf49bSBarry Smith } 350b97cf49bSBarry Smith 35140244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 352b97cf49bSBarry Smith { 353433994e6SBarry Smith PetscFunctionBegin; 3543a40ed3dSBarry Smith PetscFunctionReturn(0); 355b97cf49bSBarry Smith } 356b97cf49bSBarry Smith 35740244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 358b97cf49bSBarry Smith { 3593eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 360dfbe8321SBarry Smith PetscErrorCode ierr; 361ace3abfcSBarry Smith PetscBool flag; 362b97cf49bSBarry Smith 363433994e6SBarry Smith PetscFunctionBegin; 364b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 365d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 3660f5bd95cSBarry Smith flag = PETSC_FALSE; 367b97cf49bSBarry Smith } 368b97cf49bSBarry Smith 369b97cf49bSBarry Smith /* if the a->i are the same */ 370580bdb30SBarry Smith ierr = PetscArraycmp(a->i,b->i,A->rmap->n+1,&flag);CHKERRQ(ierr); 371b97cf49bSBarry Smith 372b97cf49bSBarry Smith /* if a->j are the same */ 373b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 374b97cf49bSBarry Smith 375b2566f29SBarry Smith ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 3763a40ed3dSBarry Smith PetscFunctionReturn(0); 377b97cf49bSBarry Smith } 378b97cf49bSBarry Smith 37940244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 3806cd91f64SBarry Smith { 381b24ad042SBarry Smith PetscInt i; 3826cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 3831a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 3846cd91f64SBarry Smith 3856cd91f64SBarry Smith PetscFunctionBegin; 386d0f46423SBarry Smith *m = A->rmap->n; 3876cd91f64SBarry Smith *ia = a->i; 3886cd91f64SBarry Smith *ja = a->j; 3896cd91f64SBarry Smith *done = PETSC_TRUE; 390d42735a1SBarry Smith if (oshift) { 391d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 392d42735a1SBarry Smith (*ja)[i]++; 393d42735a1SBarry Smith } 394d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 395d42735a1SBarry Smith } 396d42735a1SBarry Smith PetscFunctionReturn(0); 397d42735a1SBarry Smith } 398d42735a1SBarry Smith 39940244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 400d42735a1SBarry Smith { 401b24ad042SBarry Smith PetscInt i; 402d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 4031a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 404d42735a1SBarry Smith 405d42735a1SBarry Smith PetscFunctionBegin; 40670d8d27cSBarry Smith if (ia && a->i != *ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 40770d8d27cSBarry Smith if (ja && a->j != *ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 408d42735a1SBarry Smith if (oshift) { 40970d8d27cSBarry Smith if (!ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument"); 41070d8d27cSBarry Smith if (!ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument"); 411d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 412d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 413d42735a1SBarry Smith (*ja)[i]--; 414d42735a1SBarry Smith } 415d42735a1SBarry Smith } 4166cd91f64SBarry Smith PetscFunctionReturn(0); 4176cd91f64SBarry Smith } 418b97cf49bSBarry Smith 41919fd82e9SBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 42017667f90SBarry Smith { 42117667f90SBarry Smith Mat B; 42217667f90SBarry Smith PetscErrorCode ierr; 42317667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 42417667f90SBarry Smith const PetscInt *rj; 42517667f90SBarry Smith const PetscScalar *ra; 42617667f90SBarry Smith MPI_Comm comm; 42717667f90SBarry Smith 42817667f90SBarry Smith PetscFunctionBegin; 4290298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 4300298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 4310298fd71SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 43217667f90SBarry Smith 43317667f90SBarry Smith /* count the number of nonzeros per row */ 43417667f90SBarry Smith for (i=0; i<m; i++) { 4350298fd71SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 4365ee9ba1cSJed Brown for (j=0; j<len; j++) { 4375ee9ba1cSJed Brown if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 4385ee9ba1cSJed Brown } 43917667f90SBarry Smith nzeros += len; 4400f2063bfSTobin Isaac ierr = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 44117667f90SBarry Smith } 44217667f90SBarry Smith 44317667f90SBarry Smith /* malloc space for nonzeros */ 444854ce69bSBarry Smith ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr); 445854ce69bSBarry Smith ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr); 446854ce69bSBarry Smith ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr); 44717667f90SBarry Smith 44817667f90SBarry Smith nzeros = 0; 44917667f90SBarry Smith ia[0] = 0; 45017667f90SBarry Smith for (i=0; i<m; i++) { 45117667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 45217667f90SBarry Smith cnt = 0; 45317667f90SBarry Smith for (j=0; j<len; j++) { 4545ee9ba1cSJed Brown if (rj[j] != i+rstart) { /* if not diagonal */ 45517667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 45617667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 45717667f90SBarry Smith } 4585ee9ba1cSJed Brown } 45917667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 46017667f90SBarry Smith nzeros += cnt; 46117667f90SBarry Smith ia[i+1] = nzeros; 46217667f90SBarry Smith } 46317667f90SBarry Smith 46417667f90SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 46517667f90SBarry Smith ierr = MatCreate(comm,&B);CHKERRQ(ierr); 46658ec18a5SLisandro Dalcin ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 46717667f90SBarry Smith ierr = MatSetType(B,type);CHKERRQ(ierr); 46817667f90SBarry Smith ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 46917667f90SBarry Smith 470511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 47128be2f97SBarry Smith ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); 47217667f90SBarry Smith } else { 47317667f90SBarry Smith *newmat = B; 47417667f90SBarry Smith } 47517667f90SBarry Smith PetscFunctionReturn(0); 47617667f90SBarry Smith } 47717667f90SBarry Smith 478b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 4790b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 4803eda8832SBarry Smith MatGetRow_MPIAdj, 4813eda8832SBarry Smith MatRestoreRow_MPIAdj, 482b97cf49bSBarry Smith 0, 48397304618SKris Buschelman /* 4*/ 0, 484b97cf49bSBarry Smith 0, 485b97cf49bSBarry Smith 0, 486b97cf49bSBarry Smith 0, 4870b3b1487SBarry Smith 0, 4880b3b1487SBarry Smith 0, 48997304618SKris Buschelman /*10*/ 0, 4900b3b1487SBarry Smith 0, 4910b3b1487SBarry Smith 0, 4920b3b1487SBarry Smith 0, 4930b3b1487SBarry Smith 0, 49497304618SKris Buschelman /*15*/ 0, 4953eda8832SBarry Smith MatEqual_MPIAdj, 4960b3b1487SBarry Smith 0, 4970b3b1487SBarry Smith 0, 4980b3b1487SBarry Smith 0, 49997304618SKris Buschelman /*20*/ 0, 5000b3b1487SBarry Smith 0, 5013eda8832SBarry Smith MatSetOption_MPIAdj, 5020b3b1487SBarry Smith 0, 503d519adbfSMatthew Knepley /*24*/ 0, 5040b3b1487SBarry Smith 0, 5050b3b1487SBarry Smith 0, 5060b3b1487SBarry Smith 0, 5070b3b1487SBarry Smith 0, 508d519adbfSMatthew Knepley /*29*/ 0, 5090b3b1487SBarry Smith 0, 510273d9f13SBarry Smith 0, 511273d9f13SBarry Smith 0, 5120b3b1487SBarry Smith 0, 513d519adbfSMatthew Knepley /*34*/ 0, 5140b3b1487SBarry Smith 0, 5150b3b1487SBarry Smith 0, 5160b3b1487SBarry Smith 0, 5170b3b1487SBarry Smith 0, 518d519adbfSMatthew Knepley /*39*/ 0, 5197dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 5200b3b1487SBarry Smith 0, 5210b3b1487SBarry Smith 0, 5220b3b1487SBarry Smith 0, 523d519adbfSMatthew Knepley /*44*/ 0, 5240b3b1487SBarry Smith 0, 5257d68702bSBarry Smith MatShift_Basic, 5260b3b1487SBarry Smith 0, 5270b3b1487SBarry Smith 0, 528d519adbfSMatthew Knepley /*49*/ 0, 5296cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 530d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 531b97cf49bSBarry Smith 0, 532b97cf49bSBarry Smith 0, 533d519adbfSMatthew Knepley /*54*/ 0, 534b97cf49bSBarry Smith 0, 5350752156aSBarry Smith 0, 5360752156aSBarry Smith 0, 5370b3b1487SBarry Smith 0, 538d519adbfSMatthew Knepley /*59*/ 0, 539b9b97703SBarry Smith MatDestroy_MPIAdj, 540b9b97703SBarry Smith MatView_MPIAdj, 54117667f90SBarry Smith MatConvertFrom_MPIAdj, 54297304618SKris Buschelman 0, 543d519adbfSMatthew Knepley /*64*/ 0, 54497304618SKris Buschelman 0, 54597304618SKris Buschelman 0, 54697304618SKris Buschelman 0, 54797304618SKris Buschelman 0, 548d519adbfSMatthew Knepley /*69*/ 0, 54997304618SKris Buschelman 0, 55097304618SKris Buschelman 0, 55197304618SKris Buschelman 0, 55297304618SKris Buschelman 0, 553d519adbfSMatthew Knepley /*74*/ 0, 55497304618SKris Buschelman 0, 55597304618SKris Buschelman 0, 55697304618SKris Buschelman 0, 55797304618SKris Buschelman 0, 558d519adbfSMatthew Knepley /*79*/ 0, 55997304618SKris Buschelman 0, 56097304618SKris Buschelman 0, 56197304618SKris Buschelman 0, 562865e5f61SKris Buschelman 0, 563d519adbfSMatthew Knepley /*84*/ 0, 564865e5f61SKris Buschelman 0, 565865e5f61SKris Buschelman 0, 566865e5f61SKris Buschelman 0, 567865e5f61SKris Buschelman 0, 568d519adbfSMatthew Knepley /*89*/ 0, 569865e5f61SKris Buschelman 0, 570865e5f61SKris Buschelman 0, 571865e5f61SKris Buschelman 0, 572865e5f61SKris Buschelman 0, 573d519adbfSMatthew Knepley /*94*/ 0, 574865e5f61SKris Buschelman 0, 575865e5f61SKris Buschelman 0, 5763964eb88SJed Brown 0, 5773964eb88SJed Brown 0, 5783964eb88SJed Brown /*99*/ 0, 5793964eb88SJed Brown 0, 5803964eb88SJed Brown 0, 5813964eb88SJed Brown 0, 5823964eb88SJed Brown 0, 5833964eb88SJed Brown /*104*/ 0, 5843964eb88SJed Brown 0, 5853964eb88SJed Brown 0, 5863964eb88SJed Brown 0, 5873964eb88SJed Brown 0, 5883964eb88SJed Brown /*109*/ 0, 5893964eb88SJed Brown 0, 5903964eb88SJed Brown 0, 5913964eb88SJed Brown 0, 5923964eb88SJed Brown 0, 5933964eb88SJed Brown /*114*/ 0, 5943964eb88SJed Brown 0, 5953964eb88SJed Brown 0, 5963964eb88SJed Brown 0, 5973964eb88SJed Brown 0, 5983964eb88SJed Brown /*119*/ 0, 5993964eb88SJed Brown 0, 6003964eb88SJed Brown 0, 6013964eb88SJed Brown 0, 6023964eb88SJed Brown 0, 6033964eb88SJed Brown /*124*/ 0, 6043964eb88SJed Brown 0, 6053964eb88SJed Brown 0, 6063964eb88SJed Brown 0, 6077dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 6083964eb88SJed Brown /*129*/ 0, 6093964eb88SJed Brown 0, 6103964eb88SJed Brown 0, 6113964eb88SJed Brown 0, 6123964eb88SJed Brown 0, 6133964eb88SJed Brown /*134*/ 0, 6143964eb88SJed Brown 0, 6153964eb88SJed Brown 0, 6163964eb88SJed Brown 0, 6173964eb88SJed Brown 0, 6183964eb88SJed Brown /*139*/ 0, 619f9426fe0SMark Adams 0, 6203964eb88SJed Brown 0 6213964eb88SJed Brown }; 622b97cf49bSBarry Smith 623f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 624a23d5eceSKris Buschelman { 625a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 626e895ccc0SFande Kong PetscBool useedgeweights; 627dfbe8321SBarry Smith PetscErrorCode ierr; 6282515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 629b24ad042SBarry Smith PetscInt ii; 630a23d5eceSKris Buschelman #endif 631a23d5eceSKris Buschelman 632a23d5eceSKris Buschelman PetscFunctionBegin; 63326283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 63426283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 635e895ccc0SFande Kong if (values) useedgeweights = PETSC_TRUE; else useedgeweights = PETSC_FALSE; 636e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */ 637e895ccc0SFande Kong ierr = MPIU_Allreduce((int*)&useedgeweights,(int*)&b->useedgeweights,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)B));CHKERRQ(ierr); 63858ec18a5SLisandro Dalcin 6392515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 640e32f2f54SBarry 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]); 641d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 642e02043d6SBarry 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]); 643a23d5eceSKris Buschelman } 644d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 645e02043d6SBarry 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]); 646a23d5eceSKris Buschelman } 647a23d5eceSKris Buschelman #endif 64858ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 649a23d5eceSKris Buschelman 650a23d5eceSKris Buschelman b->j = j; 651a23d5eceSKris Buschelman b->i = i; 652a23d5eceSKris Buschelman b->values = values; 653a23d5eceSKris Buschelman 654d0f46423SBarry Smith b->nz = i[B->rmap->n]; 655a23d5eceSKris Buschelman b->diag = 0; 656a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 657a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 658a23d5eceSKris Buschelman 659a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 660a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 661a23d5eceSKris Buschelman PetscFunctionReturn(0); 662a23d5eceSKris Buschelman } 663b97cf49bSBarry Smith 664f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 6659a3665c6SJed Brown { 6669a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 6679a3665c6SJed Brown PetscErrorCode ierr; 6689a3665c6SJed Brown const PetscInt *ranges; 6699a3665c6SJed Brown MPI_Comm acomm,bcomm; 6709a3665c6SJed Brown MPI_Group agroup,bgroup; 6719a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 6729a3665c6SJed Brown 6739a3665c6SJed Brown PetscFunctionBegin; 6740298fd71SBarry Smith *B = NULL; 675ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr); 6769a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&size);CHKERRQ(ierr); 6779a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr); 6789a3665c6SJed Brown ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 6799a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6809a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 6819a3665c6SJed Brown } 6829a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 6839a3665c6SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 6849a3665c6SJed Brown *B = A; 6859a3665c6SJed Brown PetscFunctionReturn(0); 6869a3665c6SJed Brown } 6879a3665c6SJed Brown 688785e854fSJed Brown ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr); 6899a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6909a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 6919a3665c6SJed Brown } 6929a3665c6SJed Brown ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr); 6939a3665c6SJed Brown ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr); 6949a3665c6SJed Brown ierr = PetscFree(ranks);CHKERRQ(ierr); 6959a3665c6SJed Brown ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr); 6969a3665c6SJed Brown ierr = MPI_Group_free(&agroup);CHKERRQ(ierr); 6979a3665c6SJed Brown ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr); 6989a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 6999a3665c6SJed Brown PetscInt m,N; 7009a3665c6SJed Brown Mat_MPIAdj *b; 7010298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 7020298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 7039a3665c6SJed Brown ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr); 7049a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 7059a3665c6SJed Brown b->freeaij = PETSC_FALSE; 7069a3665c6SJed Brown ierr = MPI_Comm_free(&bcomm);CHKERRQ(ierr); 7079a3665c6SJed Brown } 7089a3665c6SJed Brown PetscFunctionReturn(0); 7099a3665c6SJed Brown } 7109a3665c6SJed Brown 711b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B) 712b44e7856SBarry Smith { 713b44e7856SBarry Smith PetscErrorCode ierr; 714b44e7856SBarry Smith PetscInt M,N,*II,*J,NZ,nz,m,nzstart,i; 715b44e7856SBarry Smith PetscInt *Values = NULL; 716b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 717b44e7856SBarry Smith PetscMPIInt mnz,mm,*allnz,*allm,size,*dispnz,*dispm; 718b44e7856SBarry Smith 719b44e7856SBarry Smith PetscFunctionBegin; 720b44e7856SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 721b44e7856SBarry Smith ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 722b44e7856SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 723b44e7856SBarry Smith nz = adj->nz; 724d4e69552SBarry Smith if (adj->i[m] != nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nz %D not correct i[m] %d",nz,adj->i[m]); 725b44e7856SBarry Smith ierr = MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 726d4e69552SBarry Smith 727b44e7856SBarry Smith ierr = PetscMPIIntCast(nz,&mnz);CHKERRQ(ierr); 728*071fcb05SBarry Smith ierr = PetscMalloc2(size,&allnz,size,&dispnz);CHKERRQ(ierr); 729d4e69552SBarry Smith ierr = MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 730b44e7856SBarry Smith dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1]; 731b44e7856SBarry Smith if (adj->values) { 732b44e7856SBarry Smith ierr = PetscMalloc1(NZ,&Values);CHKERRQ(ierr); 733b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 734b44e7856SBarry Smith } 735b44e7856SBarry Smith ierr = PetscMalloc1(NZ,&J);CHKERRQ(ierr); 736b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 737b44e7856SBarry Smith ierr = PetscFree2(allnz,dispnz);CHKERRQ(ierr); 738b44e7856SBarry Smith ierr = MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 739b44e7856SBarry Smith nzstart -= nz; 740b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 741b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] += nzstart; 742b44e7856SBarry Smith ierr = PetscMalloc1(M+1,&II);CHKERRQ(ierr); 743b44e7856SBarry Smith ierr = PetscMPIIntCast(m,&mm);CHKERRQ(ierr); 744b44e7856SBarry Smith ierr = PetscMalloc2(size,&allm,size,&dispm);CHKERRQ(ierr); 745d4e69552SBarry Smith ierr = MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 746b44e7856SBarry Smith dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1]; 747b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 748b44e7856SBarry Smith ierr = PetscFree2(allm,dispm);CHKERRQ(ierr); 749b44e7856SBarry Smith II[M] = NZ; 750b44e7856SBarry Smith /* shift the i[] values back */ 751b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] -= nzstart; 752b44e7856SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B);CHKERRQ(ierr); 753b44e7856SBarry Smith PetscFunctionReturn(0); 754b44e7856SBarry Smith } 755b44e7856SBarry Smith 7569a3665c6SJed Brown /*@ 7579a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 7589a3665c6SJed Brown 7599a3665c6SJed Brown Collective 7609a3665c6SJed Brown 7619a3665c6SJed Brown Input Arguments: 7629a3665c6SJed Brown . A - original MPIAdj matrix 7639a3665c6SJed Brown 7649a3665c6SJed Brown Output Arguments: 7650298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 7669a3665c6SJed Brown 7679a3665c6SJed Brown Level: developer 7689a3665c6SJed Brown 7699a3665c6SJed Brown Note: 7709a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 7719a3665c6SJed Brown 7729a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 7739a3665c6SJed Brown 7749a3665c6SJed Brown .seealso: MatCreateMPIAdj() 7759a3665c6SJed Brown @*/ 7769a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 7779a3665c6SJed Brown { 7789a3665c6SJed Brown PetscErrorCode ierr; 7799a3665c6SJed Brown 7809a3665c6SJed Brown PetscFunctionBegin; 7819a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 7829a3665c6SJed Brown ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 7839a3665c6SJed Brown PetscFunctionReturn(0); 7849a3665c6SJed Brown } 7859a3665c6SJed Brown 7860bad9183SKris Buschelman /*MC 787fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 7880bad9183SKris Buschelman intended for use constructing orderings and partitionings. 7890bad9183SKris Buschelman 7900bad9183SKris Buschelman Level: beginner 7910bad9183SKris Buschelman 7920bad9183SKris Buschelman .seealso: MatCreateMPIAdj 7930bad9183SKris Buschelman M*/ 7940bad9183SKris Buschelman 7958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 796273d9f13SBarry Smith { 797273d9f13SBarry Smith Mat_MPIAdj *b; 7986849ba73SBarry Smith PetscErrorCode ierr; 799273d9f13SBarry Smith 800273d9f13SBarry Smith PetscFunctionBegin; 801b00a9115SJed Brown ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 802b0a32e0cSBarry Smith B->data = (void*)b; 803273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 804273d9f13SBarry Smith B->assembled = PETSC_FALSE; 805273d9f13SBarry Smith 806bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 807bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr); 808b44e7856SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj);CHKERRQ(ierr); 80917667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 810273d9f13SBarry Smith PetscFunctionReturn(0); 811273d9f13SBarry Smith } 812273d9f13SBarry Smith 813a23d5eceSKris Buschelman /*@C 814b44e7856SBarry Smith MatMPIAdjToSeq - Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential preconditioners) 815b44e7856SBarry Smith 816d083f849SBarry Smith Logically Collective 817b44e7856SBarry Smith 818b44e7856SBarry Smith Input Parameter: 819b44e7856SBarry Smith . A - the matrix 820b44e7856SBarry Smith 821b44e7856SBarry Smith Output Parameter: 822b44e7856SBarry Smith . B - the same matrix on all processes 823b44e7856SBarry Smith 824b44e7856SBarry Smith Level: intermediate 825b44e7856SBarry Smith 826b44e7856SBarry Smith .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 827b44e7856SBarry Smith @*/ 828b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq(Mat A,Mat *B) 829b44e7856SBarry Smith { 830b44e7856SBarry Smith PetscErrorCode ierr; 831b44e7856SBarry Smith 832b44e7856SBarry Smith PetscFunctionBegin; 833b44e7856SBarry Smith ierr = PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 834b44e7856SBarry Smith PetscFunctionReturn(0); 835b44e7856SBarry Smith } 836b44e7856SBarry Smith 837b44e7856SBarry Smith /*@C 838a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 839a23d5eceSKris Buschelman 840d083f849SBarry Smith Logically Collective 841a23d5eceSKris Buschelman 842a23d5eceSKris Buschelman Input Parameters: 843a23d5eceSKris Buschelman + A - the matrix 844a23d5eceSKris Buschelman . i - the indices into j for the start of each row 845a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 846a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 847a23d5eceSKris Buschelman - values - [optional] edge weights 848a23d5eceSKris Buschelman 849a23d5eceSKris Buschelman Level: intermediate 850a23d5eceSKris Buschelman 851a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 852a23d5eceSKris Buschelman @*/ 8537087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 854273d9f13SBarry Smith { 8554ac538c5SBarry Smith PetscErrorCode ierr; 856273d9f13SBarry Smith 857273d9f13SBarry Smith PetscFunctionBegin; 8584ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 859273d9f13SBarry Smith PetscFunctionReturn(0); 860273d9f13SBarry Smith } 861273d9f13SBarry Smith 862b97cf49bSBarry Smith /*@C 8633eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 864b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 865b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 866b97cf49bSBarry Smith 867d083f849SBarry Smith Collective 868ef5ee4d1SLois Curfman McInnes 869b97cf49bSBarry Smith Input Parameters: 870c2e958feSBarry Smith + comm - MPI communicator 8710752156aSBarry Smith . m - number of local rows 87258ec18a5SLisandro Dalcin . N - number of global columns 873b97cf49bSBarry Smith . i - the indices into j for the start of each row 874329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 875ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 876329f5518SBarry Smith - values -[optional] edge weights 877b97cf49bSBarry Smith 878b97cf49bSBarry Smith Output Parameter: 879b97cf49bSBarry Smith . A - the matrix 880b97cf49bSBarry Smith 88115091d37SBarry Smith Level: intermediate 88215091d37SBarry Smith 88395452b02SPatrick Sanan Notes: 88495452b02SPatrick Sanan This matrix object does not support most matrix operations, include 8854bc6d8c0SBarry Smith MatSetValues(). 886c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 887ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 8881198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 889327e1a73SBarry Smith Should not include the matrix diagonals. 890b97cf49bSBarry Smith 891e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 892ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 893ededeb1aSvictorle 894ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 895b97cf49bSBarry Smith 896e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 897b97cf49bSBarry Smith @*/ 8987087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 899b97cf49bSBarry Smith { 900dfbe8321SBarry Smith PetscErrorCode ierr; 901b97cf49bSBarry Smith 902433994e6SBarry Smith PetscFunctionBegin; 903f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 90458ec18a5SLisandro Dalcin ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 905273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 906273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 9073a40ed3dSBarry Smith PetscFunctionReturn(0); 908b97cf49bSBarry Smith } 909b97cf49bSBarry Smith 910c2e958feSBarry Smith 911433994e6SBarry Smith 912433994e6SBarry Smith 913433994e6SBarry Smith 914433994e6SBarry Smith 915433994e6SBarry Smith 916