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; 1540244768SBarry Smith PetscInt *ncols_recv_offsets,loc,rnclos,*sadjncy,*sxadj,*svalues,isvalue; 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); 3040244768SBarry Smith ierr = PetscCalloc1(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 ???*/ 7740244768SBarry Smith if (a->values){isvalue=1;} else {isvalue=0;} 7840244768SBarry Smith /*involve a global communication */ 7940244768SBarry Smith /*ierr = MPIU_Allreduce(&isvalue,&isvalue,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);*/ 8040244768SBarry Smith if (isvalue){ierr = PetscCalloc1(Ncols_recv,&values_recv);CHKERRQ(ierr);} 8140244768SBarry Smith nroots = Ncols_send; 8240244768SBarry Smith ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 83390e1bf2SBarry Smith ierr = PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 8440244768SBarry Smith ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 8540244768SBarry Smith ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 8640244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 8740244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 8840244768SBarry Smith if (isvalue){ 8940244768SBarry Smith ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 9040244768SBarry Smith ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 9140244768SBarry Smith } 9240244768SBarry Smith ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 9340244768SBarry Smith ierr = MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr); 9440244768SBarry Smith ierr = ISGetLocalSize(icols,&icols_n);CHKERRQ(ierr); 9540244768SBarry Smith ierr = ISGetIndices(icols,&icols_indices);CHKERRQ(ierr); 9640244768SBarry Smith rnclos = 0; 9740244768SBarry Smith for (i=0; i<nlrows_is; i++){ 9840244768SBarry Smith for (j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 9940244768SBarry Smith ierr = PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc);CHKERRQ(ierr); 10040244768SBarry Smith if (loc<0){ 10140244768SBarry Smith adjncy_recv[j] = -1; 10240244768SBarry Smith if (isvalue) values_recv[j] = -1; 10340244768SBarry Smith ncols_recv[i]--; 10440244768SBarry Smith } else { 10540244768SBarry Smith rnclos++; 10640244768SBarry Smith } 10740244768SBarry Smith } 10840244768SBarry Smith } 10940244768SBarry Smith ierr = ISRestoreIndices(icols,&icols_indices);CHKERRQ(ierr); 11040244768SBarry Smith ierr = PetscCalloc1(rnclos,&sadjncy);CHKERRQ(ierr); 11140244768SBarry Smith if (isvalue) {ierr = PetscCalloc1(rnclos,&svalues);CHKERRQ(ierr);} 11240244768SBarry Smith ierr = PetscCalloc1(nlrows_is+1,&sxadj);CHKERRQ(ierr); 11340244768SBarry Smith rnclos = 0; 11440244768SBarry Smith for(i=0; i<nlrows_is; i++){ 11540244768SBarry Smith for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 11640244768SBarry Smith if (adjncy_recv[j]<0) continue; 11740244768SBarry Smith sadjncy[rnclos] = adjncy_recv[j]; 11840244768SBarry Smith if (isvalue) svalues[rnclos] = values_recv[j]; 11940244768SBarry Smith rnclos++; 12040244768SBarry Smith } 12140244768SBarry Smith } 12240244768SBarry Smith for(i=0; i<nlrows_is; i++){ 12340244768SBarry Smith sxadj[i+1] = sxadj[i]+ncols_recv[i]; 12440244768SBarry Smith } 12540244768SBarry Smith if (sadj_xadj) { *sadj_xadj = sxadj;} else { ierr = PetscFree(sxadj);CHKERRQ(ierr);} 12640244768SBarry Smith if (sadj_adjncy){ *sadj_adjncy = sadjncy;} else { ierr = PetscFree(sadjncy);CHKERRQ(ierr);} 12740244768SBarry Smith if (sadj_values){ 12840244768SBarry Smith if (isvalue) *sadj_values = svalues; else *sadj_values=0; 12940244768SBarry Smith } else { 13040244768SBarry Smith if (isvalue) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 13140244768SBarry Smith } 13240244768SBarry Smith ierr = PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv);CHKERRQ(ierr); 13340244768SBarry Smith ierr = PetscFree(adjncy_recv);CHKERRQ(ierr); 13440244768SBarry Smith if (isvalue) {ierr = PetscFree(values_recv);CHKERRQ(ierr);} 13540244768SBarry Smith PetscFunctionReturn(0); 13640244768SBarry Smith } 13740244768SBarry Smith 1387dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj_Private(Mat mat,PetscInt n,const IS irow[],const IS icol[],PetscBool subcomm,MatReuse scall,Mat *submat[]) 13940244768SBarry Smith { 14040244768SBarry Smith PetscInt i,irow_n,icol_n,*sxadj,*sadjncy,*svalues; 14140244768SBarry Smith PetscInt *indices,nindx,j,k,loc; 14240244768SBarry Smith PetscMPIInt issame; 14340244768SBarry Smith const PetscInt *irow_indices,*icol_indices; 14440244768SBarry Smith MPI_Comm scomm_row,scomm_col,scomm_mat; 14540244768SBarry Smith PetscErrorCode ierr; 14640244768SBarry Smith 14740244768SBarry Smith PetscFunctionBegin; 14840244768SBarry Smith nindx = 0; 14940244768SBarry Smith /* 15040244768SBarry Smith * Estimate a maximum number for allocating memory 15140244768SBarry Smith */ 15240244768SBarry Smith for(i=0; i<n; i++){ 15340244768SBarry Smith ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr); 15440244768SBarry Smith ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr); 15540244768SBarry Smith nindx = nindx>(irow_n+icol_n)? nindx:(irow_n+icol_n); 15640244768SBarry Smith } 15740244768SBarry Smith ierr = PetscCalloc1(nindx,&indices);CHKERRQ(ierr); 15840244768SBarry Smith /* construct a submat */ 15940244768SBarry Smith for(i=0; i<n; i++){ 16040244768SBarry Smith if (subcomm){ 16140244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)irow[i],&scomm_row);CHKERRQ(ierr); 16240244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)icol[i],&scomm_col);CHKERRQ(ierr); 16340244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,scomm_col,&issame);CHKERRQ(ierr); 16440244768SBarry 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"); 16540244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame);CHKERRQ(ierr); 16640244768SBarry 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"); 16740244768SBarry Smith } else { 16840244768SBarry Smith scomm_row = PETSC_COMM_SELF; 16940244768SBarry Smith } 17040244768SBarry Smith /*get sub-matrix data*/ 17140244768SBarry Smith sxadj=0; sadjncy=0; svalues=0; 1727dae84e0SHong Zhang ierr = MatCreateSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues);CHKERRQ(ierr); 17340244768SBarry Smith ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr); 17440244768SBarry Smith ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr); 17540244768SBarry Smith ierr = ISGetIndices(irow[i],&irow_indices);CHKERRQ(ierr); 17640244768SBarry Smith ierr = PetscMemcpy(indices,irow_indices,sizeof(PetscInt)*irow_n);CHKERRQ(ierr); 17740244768SBarry Smith ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr); 17840244768SBarry Smith ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr); 17940244768SBarry Smith ierr = PetscMemcpy(indices+irow_n,icol_indices,sizeof(PetscInt)*icol_n);CHKERRQ(ierr); 18040244768SBarry Smith ierr = ISRestoreIndices(icol[i],&icol_indices);CHKERRQ(ierr); 18140244768SBarry Smith nindx = irow_n+icol_n; 18240244768SBarry Smith ierr = PetscSortRemoveDupsInt(&nindx,indices);CHKERRQ(ierr); 18340244768SBarry Smith /* renumber columns */ 18440244768SBarry Smith for(j=0; j<irow_n; j++){ 18540244768SBarry Smith for(k=sxadj[j]; k<sxadj[j+1]; k++){ 18640244768SBarry Smith ierr = PetscFindInt(sadjncy[k],nindx,indices,&loc);CHKERRQ(ierr); 18740244768SBarry Smith #if PETSC_USE_DEBUG 18840244768SBarry Smith if (loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %d \n",sadjncy[k]); 18940244768SBarry Smith #endif 19040244768SBarry Smith sadjncy[k] = loc; 19140244768SBarry Smith } 19240244768SBarry Smith } 19340244768SBarry Smith if (scall==MAT_INITIAL_MATRIX){ 19440244768SBarry Smith ierr = MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr); 19540244768SBarry Smith } else { 19640244768SBarry Smith Mat sadj = *(submat[i]); 19740244768SBarry Smith Mat_MPIAdj *sa = (Mat_MPIAdj*)((sadj)->data); 19840244768SBarry Smith ierr = PetscObjectGetComm((PetscObject)sadj,&scomm_mat);CHKERRQ(ierr); 19940244768SBarry Smith ierr = MPI_Comm_compare(scomm_row,scomm_mat,&issame);CHKERRQ(ierr); 20040244768SBarry 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"); 20140244768SBarry Smith ierr = PetscMemcpy(sa->i,sxadj,sizeof(PetscInt)*(irow_n+1));CHKERRQ(ierr); 20240244768SBarry Smith ierr = PetscMemcpy(sa->j,sadjncy,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr); 20340244768SBarry Smith if (svalues){ierr = PetscMemcpy(sa->values,svalues,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);} 20440244768SBarry Smith ierr = PetscFree(sxadj);CHKERRQ(ierr); 20540244768SBarry Smith ierr = PetscFree(sadjncy);CHKERRQ(ierr); 20640244768SBarry Smith if (svalues) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 20740244768SBarry Smith } 20840244768SBarry Smith } 20940244768SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 21040244768SBarry Smith PetscFunctionReturn(0); 21140244768SBarry Smith } 21240244768SBarry Smith 2137dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatricesMPI_MPIAdj(Mat mat,PetscInt n, const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 21440244768SBarry Smith { 21540244768SBarry Smith PetscErrorCode ierr; 21640244768SBarry Smith /*get sub-matrices across a sub communicator */ 21740244768SBarry Smith PetscFunctionBegin; 2187dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat);CHKERRQ(ierr); 21940244768SBarry Smith PetscFunctionReturn(0); 22040244768SBarry Smith } 22140244768SBarry Smith 2227dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 22340244768SBarry Smith { 22440244768SBarry Smith PetscErrorCode ierr; 22540244768SBarry Smith 22640244768SBarry Smith PetscFunctionBegin; 22740244768SBarry Smith /*get sub-matrices based on PETSC_COMM_SELF */ 2287dae84e0SHong Zhang ierr = MatCreateSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat);CHKERRQ(ierr); 22940244768SBarry Smith PetscFunctionReturn(0); 23040244768SBarry Smith } 23140244768SBarry Smith 23240244768SBarry Smith static PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 233b97cf49bSBarry Smith { 2343eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 235dfbe8321SBarry Smith PetscErrorCode ierr; 236d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 2372dcb1b2aSMatthew Knepley const char *name; 238ce1411ecSBarry Smith PetscViewerFormat format; 239b97cf49bSBarry Smith 240433994e6SBarry Smith PetscFunctionBegin; 2413a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 242b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 243fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 2443a40ed3dSBarry Smith PetscFunctionReturn(0); 2456c4ed002SBarry Smith } else if (format == PETSC_VIEWER_ASCII_MATLAB) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported"); 2466c4ed002SBarry Smith else { 247d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 248c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 249b97cf49bSBarry Smith for (i=0; i<m; i++) { 250d0f46423SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr); 251b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 25277431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 2530752156aSBarry Smith } 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) 283d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",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); 298dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);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: 319290bbb0aSBarry Smith ierr = PetscInfo1(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; 332e32f2f54SBarry Smith if (row < 0 || row >= A->rmap->n) SETERRQ(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 */ 369d0f46423SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&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 374b2566f29SBarry Smith ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(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; 40570d8d27cSBarry Smith if (ia && a->i != *ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 40670d8d27cSBarry Smith if (ja && a->j != *ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 407d42735a1SBarry Smith if (oshift) { 40870d8d27cSBarry Smith if (!ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument"); 40970d8d27cSBarry Smith if (!ja) SETERRQ(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 /* -------------------------------------------------------------------*/ 4780b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 4793eda8832SBarry Smith MatGetRow_MPIAdj, 4803eda8832SBarry Smith MatRestoreRow_MPIAdj, 481b97cf49bSBarry Smith 0, 48297304618SKris Buschelman /* 4*/ 0, 483b97cf49bSBarry Smith 0, 484b97cf49bSBarry Smith 0, 485b97cf49bSBarry Smith 0, 4860b3b1487SBarry Smith 0, 4870b3b1487SBarry Smith 0, 48897304618SKris Buschelman /*10*/ 0, 4890b3b1487SBarry Smith 0, 4900b3b1487SBarry Smith 0, 4910b3b1487SBarry Smith 0, 4920b3b1487SBarry Smith 0, 49397304618SKris Buschelman /*15*/ 0, 4943eda8832SBarry Smith MatEqual_MPIAdj, 4950b3b1487SBarry Smith 0, 4960b3b1487SBarry Smith 0, 4970b3b1487SBarry Smith 0, 49897304618SKris Buschelman /*20*/ 0, 4990b3b1487SBarry Smith 0, 5003eda8832SBarry Smith MatSetOption_MPIAdj, 5010b3b1487SBarry Smith 0, 502d519adbfSMatthew Knepley /*24*/ 0, 5030b3b1487SBarry Smith 0, 5040b3b1487SBarry Smith 0, 5050b3b1487SBarry Smith 0, 5060b3b1487SBarry Smith 0, 507d519adbfSMatthew Knepley /*29*/ 0, 5080b3b1487SBarry Smith 0, 509273d9f13SBarry Smith 0, 510273d9f13SBarry Smith 0, 5110b3b1487SBarry Smith 0, 512d519adbfSMatthew Knepley /*34*/ 0, 5130b3b1487SBarry Smith 0, 5140b3b1487SBarry Smith 0, 5150b3b1487SBarry Smith 0, 5160b3b1487SBarry Smith 0, 517d519adbfSMatthew Knepley /*39*/ 0, 5187dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 5190b3b1487SBarry Smith 0, 5200b3b1487SBarry Smith 0, 5210b3b1487SBarry Smith 0, 522d519adbfSMatthew Knepley /*44*/ 0, 5230b3b1487SBarry Smith 0, 5247d68702bSBarry Smith MatShift_Basic, 5250b3b1487SBarry Smith 0, 5260b3b1487SBarry Smith 0, 527d519adbfSMatthew Knepley /*49*/ 0, 5286cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 529d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 530b97cf49bSBarry Smith 0, 531b97cf49bSBarry Smith 0, 532d519adbfSMatthew Knepley /*54*/ 0, 533b97cf49bSBarry Smith 0, 5340752156aSBarry Smith 0, 5350752156aSBarry Smith 0, 5360b3b1487SBarry Smith 0, 537d519adbfSMatthew Knepley /*59*/ 0, 538b9b97703SBarry Smith MatDestroy_MPIAdj, 539b9b97703SBarry Smith MatView_MPIAdj, 54017667f90SBarry Smith MatConvertFrom_MPIAdj, 54197304618SKris Buschelman 0, 542d519adbfSMatthew Knepley /*64*/ 0, 54397304618SKris Buschelman 0, 54497304618SKris Buschelman 0, 54597304618SKris Buschelman 0, 54697304618SKris Buschelman 0, 547d519adbfSMatthew Knepley /*69*/ 0, 54897304618SKris Buschelman 0, 54997304618SKris Buschelman 0, 55097304618SKris Buschelman 0, 55197304618SKris Buschelman 0, 552d519adbfSMatthew Knepley /*74*/ 0, 55397304618SKris Buschelman 0, 55497304618SKris Buschelman 0, 55597304618SKris Buschelman 0, 55697304618SKris Buschelman 0, 557d519adbfSMatthew Knepley /*79*/ 0, 55897304618SKris Buschelman 0, 55997304618SKris Buschelman 0, 56097304618SKris Buschelman 0, 561865e5f61SKris Buschelman 0, 562d519adbfSMatthew Knepley /*84*/ 0, 563865e5f61SKris Buschelman 0, 564865e5f61SKris Buschelman 0, 565865e5f61SKris Buschelman 0, 566865e5f61SKris Buschelman 0, 567d519adbfSMatthew Knepley /*89*/ 0, 568865e5f61SKris Buschelman 0, 569865e5f61SKris Buschelman 0, 570865e5f61SKris Buschelman 0, 571865e5f61SKris Buschelman 0, 572d519adbfSMatthew Knepley /*94*/ 0, 573865e5f61SKris Buschelman 0, 574865e5f61SKris Buschelman 0, 5753964eb88SJed Brown 0, 5763964eb88SJed Brown 0, 5773964eb88SJed Brown /*99*/ 0, 5783964eb88SJed Brown 0, 5793964eb88SJed Brown 0, 5803964eb88SJed Brown 0, 5813964eb88SJed Brown 0, 5823964eb88SJed Brown /*104*/ 0, 5833964eb88SJed Brown 0, 5843964eb88SJed Brown 0, 5853964eb88SJed Brown 0, 5863964eb88SJed Brown 0, 5873964eb88SJed Brown /*109*/ 0, 5883964eb88SJed Brown 0, 5893964eb88SJed Brown 0, 5903964eb88SJed Brown 0, 5913964eb88SJed Brown 0, 5923964eb88SJed Brown /*114*/ 0, 5933964eb88SJed Brown 0, 5943964eb88SJed Brown 0, 5953964eb88SJed Brown 0, 5963964eb88SJed Brown 0, 5973964eb88SJed Brown /*119*/ 0, 5983964eb88SJed Brown 0, 5993964eb88SJed Brown 0, 6003964eb88SJed Brown 0, 6013964eb88SJed Brown 0, 6023964eb88SJed Brown /*124*/ 0, 6033964eb88SJed Brown 0, 6043964eb88SJed Brown 0, 6053964eb88SJed Brown 0, 6067dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 6073964eb88SJed Brown /*129*/ 0, 6083964eb88SJed Brown 0, 6093964eb88SJed Brown 0, 6103964eb88SJed Brown 0, 6113964eb88SJed Brown 0, 6123964eb88SJed Brown /*134*/ 0, 6133964eb88SJed Brown 0, 6143964eb88SJed Brown 0, 6153964eb88SJed Brown 0, 6163964eb88SJed Brown 0, 6173964eb88SJed Brown /*139*/ 0, 618f9426fe0SMark Adams 0, 6193964eb88SJed Brown 0 6203964eb88SJed Brown }; 621b97cf49bSBarry Smith 622f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 623a23d5eceSKris Buschelman { 624a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 625dfbe8321SBarry Smith PetscErrorCode ierr; 6262515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 627b24ad042SBarry Smith PetscInt ii; 628a23d5eceSKris Buschelman #endif 629a23d5eceSKris Buschelman 630a23d5eceSKris Buschelman PetscFunctionBegin; 63126283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 63226283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 63358ec18a5SLisandro Dalcin 6342515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 635e32f2f54SBarry 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]); 636d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 637e02043d6SBarry 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]); 638a23d5eceSKris Buschelman } 639d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 640e02043d6SBarry 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]); 641a23d5eceSKris Buschelman } 642a23d5eceSKris Buschelman #endif 64358ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 644a23d5eceSKris Buschelman 645a23d5eceSKris Buschelman b->j = j; 646a23d5eceSKris Buschelman b->i = i; 647a23d5eceSKris Buschelman b->values = values; 648a23d5eceSKris Buschelman 649d0f46423SBarry Smith b->nz = i[B->rmap->n]; 650a23d5eceSKris Buschelman b->diag = 0; 651a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 652a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 653a23d5eceSKris Buschelman 654a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 655a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 656a23d5eceSKris Buschelman PetscFunctionReturn(0); 657a23d5eceSKris Buschelman } 658b97cf49bSBarry Smith 659f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 6609a3665c6SJed Brown { 6619a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 6629a3665c6SJed Brown PetscErrorCode ierr; 6639a3665c6SJed Brown const PetscInt *ranges; 6649a3665c6SJed Brown MPI_Comm acomm,bcomm; 6659a3665c6SJed Brown MPI_Group agroup,bgroup; 6669a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 6679a3665c6SJed Brown 6689a3665c6SJed Brown PetscFunctionBegin; 6690298fd71SBarry Smith *B = NULL; 670ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr); 6719a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&size);CHKERRQ(ierr); 6729a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr); 6739a3665c6SJed Brown ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 6749a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6759a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 6769a3665c6SJed Brown } 6779a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 6789a3665c6SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 6799a3665c6SJed Brown *B = A; 6809a3665c6SJed Brown PetscFunctionReturn(0); 6819a3665c6SJed Brown } 6829a3665c6SJed Brown 683785e854fSJed Brown ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr); 6849a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6859a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 6869a3665c6SJed Brown } 6879a3665c6SJed Brown ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr); 6889a3665c6SJed Brown ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr); 6899a3665c6SJed Brown ierr = PetscFree(ranks);CHKERRQ(ierr); 6909a3665c6SJed Brown ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr); 6919a3665c6SJed Brown ierr = MPI_Group_free(&agroup);CHKERRQ(ierr); 6929a3665c6SJed Brown ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr); 6939a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 6949a3665c6SJed Brown PetscInt m,N; 6959a3665c6SJed Brown Mat_MPIAdj *b; 6960298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 6970298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 6989a3665c6SJed Brown ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr); 6999a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 7009a3665c6SJed Brown b->freeaij = PETSC_FALSE; 7019a3665c6SJed Brown ierr = MPI_Comm_free(&bcomm);CHKERRQ(ierr); 7029a3665c6SJed Brown } 7039a3665c6SJed Brown PetscFunctionReturn(0); 7049a3665c6SJed Brown } 7059a3665c6SJed Brown 706b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B) 707b44e7856SBarry Smith { 708b44e7856SBarry Smith PetscErrorCode ierr; 709b44e7856SBarry Smith PetscInt M,N,*II,*J,NZ,nz,m,nzstart,i; 710b44e7856SBarry Smith PetscInt *Values = NULL; 711b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj*)A->data; 712b44e7856SBarry Smith PetscMPIInt mnz,mm,*allnz,*allm,size,*dispnz,*dispm; 713b44e7856SBarry Smith 714b44e7856SBarry Smith PetscFunctionBegin; 715b44e7856SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 716b44e7856SBarry Smith ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 717b44e7856SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 718b44e7856SBarry Smith nz = adj->nz; 719*d4e69552SBarry Smith if (adj->i[m] != nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nz %D not correct i[m] %d",nz,adj->i[m]); 720b44e7856SBarry Smith ierr = MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 721*d4e69552SBarry Smith 722b44e7856SBarry Smith ierr = PetscMPIIntCast(nz,&mnz);CHKERRQ(ierr); 723*d4e69552SBarry Smith ierr = PetscCalloc2(size,&allnz,size,&dispnz);CHKERRQ(ierr); 724*d4e69552SBarry Smith ierr = MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 725b44e7856SBarry Smith dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1]; 726b44e7856SBarry Smith if (adj->values) { 727b44e7856SBarry Smith ierr = PetscMalloc1(NZ,&Values);CHKERRQ(ierr); 728b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 729b44e7856SBarry Smith } 730b44e7856SBarry Smith ierr = PetscMalloc1(NZ,&J);CHKERRQ(ierr); 731b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 732b44e7856SBarry Smith ierr = PetscFree2(allnz,dispnz);CHKERRQ(ierr); 733b44e7856SBarry Smith ierr = MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 734b44e7856SBarry Smith nzstart -= nz; 735b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 736b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] += nzstart; 737b44e7856SBarry Smith ierr = PetscMalloc1(M+1,&II);CHKERRQ(ierr); 738b44e7856SBarry Smith ierr = PetscMPIIntCast(m,&mm);CHKERRQ(ierr); 739b44e7856SBarry Smith ierr = PetscMalloc2(size,&allm,size,&dispm);CHKERRQ(ierr); 740*d4e69552SBarry Smith ierr = MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 741b44e7856SBarry Smith dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1]; 742b44e7856SBarry Smith ierr = MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 743b44e7856SBarry Smith ierr = PetscFree2(allm,dispm);CHKERRQ(ierr); 744b44e7856SBarry Smith II[M] = NZ; 745b44e7856SBarry Smith /* shift the i[] values back */ 746b44e7856SBarry Smith for (i=0; i<m; i++) adj->i[i] -= nzstart; 747b44e7856SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B);CHKERRQ(ierr); 748b44e7856SBarry Smith PetscFunctionReturn(0); 749b44e7856SBarry Smith } 750b44e7856SBarry Smith 7519a3665c6SJed Brown /*@ 7529a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 7539a3665c6SJed Brown 7549a3665c6SJed Brown Collective 7559a3665c6SJed Brown 7569a3665c6SJed Brown Input Arguments: 7579a3665c6SJed Brown . A - original MPIAdj matrix 7589a3665c6SJed Brown 7599a3665c6SJed Brown Output Arguments: 7600298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 7619a3665c6SJed Brown 7629a3665c6SJed Brown Level: developer 7639a3665c6SJed Brown 7649a3665c6SJed Brown Note: 7659a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 7669a3665c6SJed Brown 7679a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 7689a3665c6SJed Brown 7699a3665c6SJed Brown .seealso: MatCreateMPIAdj() 7709a3665c6SJed Brown @*/ 7719a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 7729a3665c6SJed Brown { 7739a3665c6SJed Brown PetscErrorCode ierr; 7749a3665c6SJed Brown 7759a3665c6SJed Brown PetscFunctionBegin; 7769a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 7779a3665c6SJed Brown ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 7789a3665c6SJed Brown PetscFunctionReturn(0); 7799a3665c6SJed Brown } 7809a3665c6SJed Brown 7810bad9183SKris Buschelman /*MC 782fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 7830bad9183SKris Buschelman intended for use constructing orderings and partitionings. 7840bad9183SKris Buschelman 7850bad9183SKris Buschelman Level: beginner 7860bad9183SKris Buschelman 7870bad9183SKris Buschelman .seealso: MatCreateMPIAdj 7880bad9183SKris Buschelman M*/ 7890bad9183SKris Buschelman 7908cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 791273d9f13SBarry Smith { 792273d9f13SBarry Smith Mat_MPIAdj *b; 7936849ba73SBarry Smith PetscErrorCode ierr; 794273d9f13SBarry Smith 795273d9f13SBarry Smith PetscFunctionBegin; 796b00a9115SJed Brown ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 797b0a32e0cSBarry Smith B->data = (void*)b; 798273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 799273d9f13SBarry Smith B->assembled = PETSC_FALSE; 800273d9f13SBarry Smith 801bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 802bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr); 803b44e7856SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj);CHKERRQ(ierr); 80417667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 805273d9f13SBarry Smith PetscFunctionReturn(0); 806273d9f13SBarry Smith } 807273d9f13SBarry Smith 808a23d5eceSKris Buschelman /*@C 809b44e7856SBarry Smith MatMPIAdjToSeq- Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential preconditioners) 810b44e7856SBarry Smith 811b44e7856SBarry Smith Logically Collective on MPI_Comm 812b44e7856SBarry Smith 813b44e7856SBarry Smith Input Parameter: 814b44e7856SBarry Smith . A - the matrix 815b44e7856SBarry Smith 816b44e7856SBarry Smith Output Parameter: 817b44e7856SBarry Smith . B - the same matrix on all processes 818b44e7856SBarry Smith 819b44e7856SBarry Smith Level: intermediate 820b44e7856SBarry Smith 821b44e7856SBarry Smith .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 822b44e7856SBarry Smith @*/ 823b44e7856SBarry Smith PetscErrorCode MatMPIAdjToSeq(Mat A,Mat *B) 824b44e7856SBarry Smith { 825b44e7856SBarry Smith PetscErrorCode ierr; 826b44e7856SBarry Smith 827b44e7856SBarry Smith PetscFunctionBegin; 828b44e7856SBarry Smith ierr = PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 829b44e7856SBarry Smith PetscFunctionReturn(0); 830b44e7856SBarry Smith } 831b44e7856SBarry Smith 832b44e7856SBarry Smith /*@C 833a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 834a23d5eceSKris Buschelman 8353f9fe445SBarry Smith Logically Collective on MPI_Comm 836a23d5eceSKris Buschelman 837a23d5eceSKris Buschelman Input Parameters: 838a23d5eceSKris Buschelman + A - the matrix 839a23d5eceSKris Buschelman . i - the indices into j for the start of each row 840a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 841a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 842a23d5eceSKris Buschelman - values - [optional] edge weights 843a23d5eceSKris Buschelman 844a23d5eceSKris Buschelman Level: intermediate 845a23d5eceSKris Buschelman 846a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 847a23d5eceSKris Buschelman @*/ 8487087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 849273d9f13SBarry Smith { 8504ac538c5SBarry Smith PetscErrorCode ierr; 851273d9f13SBarry Smith 852273d9f13SBarry Smith PetscFunctionBegin; 8534ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 854273d9f13SBarry Smith PetscFunctionReturn(0); 855273d9f13SBarry Smith } 856273d9f13SBarry Smith 857b97cf49bSBarry Smith /*@C 8583eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 859b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 860b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 861b97cf49bSBarry Smith 862ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 863ef5ee4d1SLois Curfman McInnes 864b97cf49bSBarry Smith Input Parameters: 865c2e958feSBarry Smith + comm - MPI communicator 8660752156aSBarry Smith . m - number of local rows 86758ec18a5SLisandro Dalcin . N - number of global columns 868b97cf49bSBarry Smith . i - the indices into j for the start of each row 869329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 870ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 871329f5518SBarry Smith - values -[optional] edge weights 872b97cf49bSBarry Smith 873b97cf49bSBarry Smith Output Parameter: 874b97cf49bSBarry Smith . A - the matrix 875b97cf49bSBarry Smith 87615091d37SBarry Smith Level: intermediate 87715091d37SBarry Smith 8784bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 8794bc6d8c0SBarry Smith MatSetValues(). 880c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 881ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 8821198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 883327e1a73SBarry Smith Should not include the matrix diagonals. 884b97cf49bSBarry Smith 885e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 886ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 887ededeb1aSvictorle 888ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 889b97cf49bSBarry Smith 890e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 891b97cf49bSBarry Smith @*/ 8927087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 893b97cf49bSBarry Smith { 894dfbe8321SBarry Smith PetscErrorCode ierr; 895b97cf49bSBarry Smith 896433994e6SBarry Smith PetscFunctionBegin; 897f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 89858ec18a5SLisandro Dalcin ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 899273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 900273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 9013a40ed3dSBarry Smith PetscFunctionReturn(0); 902b97cf49bSBarry Smith } 903b97cf49bSBarry Smith 904c2e958feSBarry Smith 905433994e6SBarry Smith 906433994e6SBarry Smith 907433994e6SBarry Smith 908433994e6SBarry Smith 909433994e6SBarry Smith 910