xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision e895ccc02404e65fbacee93aeb4f80c54d3f2a09)
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;
15*e895ccc0SFande 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);
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 ???*/
77*e895ccc0SFande 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);
85*e895ccc0SFande 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;
99*e895ccc0SFande 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);
108*e895ccc0SFande 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];
115*e895ccc0SFande 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){
125*e895ccc0SFande Kong     if (a->useedgeweights) *sadj_values = svalues; else *sadj_values=0;
12640244768SBarry Smith   } else {
127*e895ccc0SFande 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);
131*e895ccc0SFande 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   }
15440244768SBarry Smith   ierr = PetscCalloc1(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);
17340244768SBarry Smith     ierr = PetscMemcpy(indices,irow_indices,sizeof(PetscInt)*irow_n);CHKERRQ(ierr);
17440244768SBarry Smith     ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr);
17540244768SBarry Smith     ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr);
17640244768SBarry Smith     ierr = PetscMemcpy(indices+irow_n,icol_indices,sizeof(PetscInt)*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");
19840244768SBarry Smith        ierr = PetscMemcpy(sa->i,sxadj,sizeof(PetscInt)*(irow_n+1));CHKERRQ(ierr);
19940244768SBarry Smith        ierr = PetscMemcpy(sa->j,sadjncy,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);
20040244768SBarry Smith        if (svalues){ierr = PetscMemcpy(sa->values,svalues,sizeof(PetscInt)*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++) {
24977431f27SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr);
2500752156aSBarry Smith       }
251b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
252b97cf49bSBarry Smith     }
253d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
254b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
255c5e4d11fSDmitry Karpeev     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
2567b23a99aSBarry Smith   }
2573a40ed3dSBarry Smith   PetscFunctionReturn(0);
258b97cf49bSBarry Smith }
259b97cf49bSBarry Smith 
26040244768SBarry Smith static PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer)
261b97cf49bSBarry Smith {
262dfbe8321SBarry Smith   PetscErrorCode ierr;
263ace3abfcSBarry Smith   PetscBool      iascii;
264b97cf49bSBarry Smith 
265433994e6SBarry Smith   PetscFunctionBegin;
266251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
26732077d6dSBarry Smith   if (iascii) {
2683eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
269b97cf49bSBarry Smith   }
2703a40ed3dSBarry Smith   PetscFunctionReturn(0);
271b97cf49bSBarry Smith }
272b97cf49bSBarry Smith 
27340244768SBarry Smith static PetscErrorCode MatDestroy_MPIAdj(Mat mat)
274b97cf49bSBarry Smith {
2753eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)mat->data;
276dfbe8321SBarry Smith   PetscErrorCode ierr;
27794d884c6SBarry Smith 
27894d884c6SBarry Smith   PetscFunctionBegin;
279aa482453SBarry Smith #if defined(PETSC_USE_LOG)
280d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz);
281b97cf49bSBarry Smith #endif
28205b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
2830400557aSBarry Smith   if (a->freeaij) {
28414196391SBarry Smith     if (a->freeaijwithfree) {
28514196391SBarry Smith       if (a->i) free(a->i);
28614196391SBarry Smith       if (a->j) free(a->j);
28714196391SBarry Smith     } else {
288606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
289606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
29005b42c5fSBarry Smith       ierr = PetscFree(a->values);CHKERRQ(ierr);
29114196391SBarry Smith     }
2920400557aSBarry Smith   }
2932b1d8763SJed Brown   ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
294bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
295dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
296bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr);
297bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr);
2983a40ed3dSBarry Smith   PetscFunctionReturn(0);
299b97cf49bSBarry Smith }
300b97cf49bSBarry Smith 
30140244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg)
302b97cf49bSBarry Smith {
3033eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
30463ba0a88SBarry Smith   PetscErrorCode ierr;
305b97cf49bSBarry Smith 
306433994e6SBarry Smith   PetscFunctionBegin;
30712c028f9SKris Buschelman   switch (op) {
30877e54ba9SKris Buschelman   case MAT_SYMMETRIC:
30912c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
3109a4540c5SBarry Smith   case MAT_HERMITIAN:
3114e0d8c25SBarry Smith     a->symmetric = flg;
3129a4540c5SBarry Smith     break;
3139a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
3149a4540c5SBarry Smith     break;
31512c028f9SKris Buschelman   default:
316290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
31712c028f9SKris Buschelman     break;
318b97cf49bSBarry Smith   }
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
320b97cf49bSBarry Smith }
321b97cf49bSBarry Smith 
32240244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
323b97cf49bSBarry Smith {
3243eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
3252b1d8763SJed Brown   PetscErrorCode ierr;
326b97cf49bSBarry Smith 
327433994e6SBarry Smith   PetscFunctionBegin;
328d0f46423SBarry Smith   row -= A->rmap->rstart;
329e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
330b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
3312b1d8763SJed Brown   if (v) {
3322b1d8763SJed Brown     PetscInt j;
3332b1d8763SJed Brown     if (a->rowvalues_alloc < *nz) {
3342b1d8763SJed Brown       ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
3352b1d8763SJed Brown       a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz);
336785e854fSJed Brown       ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr);
337b97cf49bSBarry Smith     }
33891f8cafeSFande Kong     for (j=0; j<*nz; j++){
33991f8cafeSFande Kong       a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0;
34091f8cafeSFande Kong     }
3412b1d8763SJed Brown     *v = (*nz) ? a->rowvalues : NULL;
342b97cf49bSBarry Smith   }
34392b6f2a0SJed Brown   if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
3443a40ed3dSBarry Smith   PetscFunctionReturn(0);
345b97cf49bSBarry Smith }
346b97cf49bSBarry Smith 
34740244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
348b97cf49bSBarry Smith {
349433994e6SBarry Smith   PetscFunctionBegin;
3503a40ed3dSBarry Smith   PetscFunctionReturn(0);
351b97cf49bSBarry Smith }
352b97cf49bSBarry Smith 
35340244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
354b97cf49bSBarry Smith {
3553eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data;
356dfbe8321SBarry Smith   PetscErrorCode ierr;
357ace3abfcSBarry Smith   PetscBool      flag;
358b97cf49bSBarry Smith 
359433994e6SBarry Smith   PetscFunctionBegin;
360b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
361d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
3620f5bd95cSBarry Smith     flag = PETSC_FALSE;
363b97cf49bSBarry Smith   }
364b97cf49bSBarry Smith 
365b97cf49bSBarry Smith   /* if the a->i are the same */
366d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
367b97cf49bSBarry Smith 
368b97cf49bSBarry Smith   /* if a->j are the same */
369b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
370b97cf49bSBarry Smith 
371b2566f29SBarry Smith   ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
3723a40ed3dSBarry Smith   PetscFunctionReturn(0);
373b97cf49bSBarry Smith }
374b97cf49bSBarry Smith 
37540244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
3766cd91f64SBarry Smith {
377b24ad042SBarry Smith   PetscInt   i;
3786cd91f64SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
3791a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
3806cd91f64SBarry Smith 
3816cd91f64SBarry Smith   PetscFunctionBegin;
382d0f46423SBarry Smith   *m    = A->rmap->n;
3836cd91f64SBarry Smith   *ia   = a->i;
3846cd91f64SBarry Smith   *ja   = a->j;
3856cd91f64SBarry Smith   *done = PETSC_TRUE;
386d42735a1SBarry Smith   if (oshift) {
387d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
388d42735a1SBarry Smith       (*ja)[i]++;
389d42735a1SBarry Smith     }
390d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
391d42735a1SBarry Smith   }
392d42735a1SBarry Smith   PetscFunctionReturn(0);
393d42735a1SBarry Smith }
394d42735a1SBarry Smith 
39540244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
396d42735a1SBarry Smith {
397b24ad042SBarry Smith   PetscInt   i;
398d42735a1SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
3991a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
400d42735a1SBarry Smith 
401d42735a1SBarry Smith   PetscFunctionBegin;
40270d8d27cSBarry Smith   if (ia && a->i != *ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
40370d8d27cSBarry Smith   if (ja && a->j != *ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
404d42735a1SBarry Smith   if (oshift) {
40570d8d27cSBarry Smith     if (!ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument");
40670d8d27cSBarry Smith     if (!ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument");
407d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
408d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
409d42735a1SBarry Smith       (*ja)[i]--;
410d42735a1SBarry Smith     }
411d42735a1SBarry Smith   }
4126cd91f64SBarry Smith   PetscFunctionReturn(0);
4136cd91f64SBarry Smith }
414b97cf49bSBarry Smith 
41519fd82e9SBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
41617667f90SBarry Smith {
41717667f90SBarry Smith   Mat               B;
41817667f90SBarry Smith   PetscErrorCode    ierr;
41917667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
42017667f90SBarry Smith   const PetscInt    *rj;
42117667f90SBarry Smith   const PetscScalar *ra;
42217667f90SBarry Smith   MPI_Comm          comm;
42317667f90SBarry Smith 
42417667f90SBarry Smith   PetscFunctionBegin;
4250298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
4260298fd71SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
4270298fd71SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
42817667f90SBarry Smith 
42917667f90SBarry Smith   /* count the number of nonzeros per row */
43017667f90SBarry Smith   for (i=0; i<m; i++) {
4310298fd71SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
4325ee9ba1cSJed Brown     for (j=0; j<len; j++) {
4335ee9ba1cSJed Brown       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
4345ee9ba1cSJed Brown     }
43517667f90SBarry Smith     nzeros += len;
4360f2063bfSTobin Isaac     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
43717667f90SBarry Smith   }
43817667f90SBarry Smith 
43917667f90SBarry Smith   /* malloc space for nonzeros */
440854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr);
441854ce69bSBarry Smith   ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr);
442854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr);
44317667f90SBarry Smith 
44417667f90SBarry Smith   nzeros = 0;
44517667f90SBarry Smith   ia[0]  = 0;
44617667f90SBarry Smith   for (i=0; i<m; i++) {
44717667f90SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
44817667f90SBarry Smith     cnt  = 0;
44917667f90SBarry Smith     for (j=0; j<len; j++) {
4505ee9ba1cSJed Brown       if (rj[j] != i+rstart) { /* if not diagonal */
45117667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
45217667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
45317667f90SBarry Smith       }
4545ee9ba1cSJed Brown     }
45517667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
45617667f90SBarry Smith     nzeros += cnt;
45717667f90SBarry Smith     ia[i+1] = nzeros;
45817667f90SBarry Smith   }
45917667f90SBarry Smith 
46017667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
46117667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
46258ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
46317667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
46417667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
46517667f90SBarry Smith 
466511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
46728be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
46817667f90SBarry Smith   } else {
46917667f90SBarry Smith     *newmat = B;
47017667f90SBarry Smith   }
47117667f90SBarry Smith   PetscFunctionReturn(0);
47217667f90SBarry Smith }
47317667f90SBarry Smith 
474b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
4750b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
4763eda8832SBarry Smith                                        MatGetRow_MPIAdj,
4773eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
478b97cf49bSBarry Smith                                        0,
47997304618SKris Buschelman                                 /* 4*/ 0,
480b97cf49bSBarry Smith                                        0,
481b97cf49bSBarry Smith                                        0,
482b97cf49bSBarry Smith                                        0,
4830b3b1487SBarry Smith                                        0,
4840b3b1487SBarry Smith                                        0,
48597304618SKris Buschelman                                 /*10*/ 0,
4860b3b1487SBarry Smith                                        0,
4870b3b1487SBarry Smith                                        0,
4880b3b1487SBarry Smith                                        0,
4890b3b1487SBarry Smith                                        0,
49097304618SKris Buschelman                                 /*15*/ 0,
4913eda8832SBarry Smith                                        MatEqual_MPIAdj,
4920b3b1487SBarry Smith                                        0,
4930b3b1487SBarry Smith                                        0,
4940b3b1487SBarry Smith                                        0,
49597304618SKris Buschelman                                 /*20*/ 0,
4960b3b1487SBarry Smith                                        0,
4973eda8832SBarry Smith                                        MatSetOption_MPIAdj,
4980b3b1487SBarry Smith                                        0,
499d519adbfSMatthew Knepley                                 /*24*/ 0,
5000b3b1487SBarry Smith                                        0,
5010b3b1487SBarry Smith                                        0,
5020b3b1487SBarry Smith                                        0,
5030b3b1487SBarry Smith                                        0,
504d519adbfSMatthew Knepley                                 /*29*/ 0,
5050b3b1487SBarry Smith                                        0,
506273d9f13SBarry Smith                                        0,
507273d9f13SBarry Smith                                        0,
5080b3b1487SBarry Smith                                        0,
509d519adbfSMatthew Knepley                                 /*34*/ 0,
5100b3b1487SBarry Smith                                        0,
5110b3b1487SBarry Smith                                        0,
5120b3b1487SBarry Smith                                        0,
5130b3b1487SBarry Smith                                        0,
514d519adbfSMatthew Knepley                                 /*39*/ 0,
5157dae84e0SHong Zhang                                        MatCreateSubMatrices_MPIAdj,
5160b3b1487SBarry Smith                                        0,
5170b3b1487SBarry Smith                                        0,
5180b3b1487SBarry Smith                                        0,
519d519adbfSMatthew Knepley                                 /*44*/ 0,
5200b3b1487SBarry Smith                                        0,
5217d68702bSBarry Smith                                        MatShift_Basic,
5220b3b1487SBarry Smith                                        0,
5230b3b1487SBarry Smith                                        0,
524d519adbfSMatthew Knepley                                 /*49*/ 0,
5256cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
526d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
527b97cf49bSBarry Smith                                        0,
528b97cf49bSBarry Smith                                        0,
529d519adbfSMatthew Knepley                                 /*54*/ 0,
530b97cf49bSBarry Smith                                        0,
5310752156aSBarry Smith                                        0,
5320752156aSBarry Smith                                        0,
5330b3b1487SBarry Smith                                        0,
534d519adbfSMatthew Knepley                                 /*59*/ 0,
535b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
536b9b97703SBarry Smith                                        MatView_MPIAdj,
53717667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
53897304618SKris Buschelman                                        0,
539d519adbfSMatthew Knepley                                 /*64*/ 0,
54097304618SKris Buschelman                                        0,
54197304618SKris Buschelman                                        0,
54297304618SKris Buschelman                                        0,
54397304618SKris Buschelman                                        0,
544d519adbfSMatthew Knepley                                 /*69*/ 0,
54597304618SKris Buschelman                                        0,
54697304618SKris Buschelman                                        0,
54797304618SKris Buschelman                                        0,
54897304618SKris Buschelman                                        0,
549d519adbfSMatthew Knepley                                 /*74*/ 0,
55097304618SKris Buschelman                                        0,
55197304618SKris Buschelman                                        0,
55297304618SKris Buschelman                                        0,
55397304618SKris Buschelman                                        0,
554d519adbfSMatthew Knepley                                 /*79*/ 0,
55597304618SKris Buschelman                                        0,
55697304618SKris Buschelman                                        0,
55797304618SKris Buschelman                                        0,
558865e5f61SKris Buschelman                                        0,
559d519adbfSMatthew Knepley                                 /*84*/ 0,
560865e5f61SKris Buschelman                                        0,
561865e5f61SKris Buschelman                                        0,
562865e5f61SKris Buschelman                                        0,
563865e5f61SKris Buschelman                                        0,
564d519adbfSMatthew Knepley                                 /*89*/ 0,
565865e5f61SKris Buschelman                                        0,
566865e5f61SKris Buschelman                                        0,
567865e5f61SKris Buschelman                                        0,
568865e5f61SKris Buschelman                                        0,
569d519adbfSMatthew Knepley                                 /*94*/ 0,
570865e5f61SKris Buschelman                                        0,
571865e5f61SKris Buschelman                                        0,
5723964eb88SJed Brown                                        0,
5733964eb88SJed Brown                                        0,
5743964eb88SJed Brown                                 /*99*/ 0,
5753964eb88SJed Brown                                        0,
5763964eb88SJed Brown                                        0,
5773964eb88SJed Brown                                        0,
5783964eb88SJed Brown                                        0,
5793964eb88SJed Brown                                /*104*/ 0,
5803964eb88SJed Brown                                        0,
5813964eb88SJed Brown                                        0,
5823964eb88SJed Brown                                        0,
5833964eb88SJed Brown                                        0,
5843964eb88SJed Brown                                /*109*/ 0,
5853964eb88SJed Brown                                        0,
5863964eb88SJed Brown                                        0,
5873964eb88SJed Brown                                        0,
5883964eb88SJed Brown                                        0,
5893964eb88SJed Brown                                /*114*/ 0,
5903964eb88SJed Brown                                        0,
5913964eb88SJed Brown                                        0,
5923964eb88SJed Brown                                        0,
5933964eb88SJed Brown                                        0,
5943964eb88SJed Brown                                /*119*/ 0,
5953964eb88SJed Brown                                        0,
5963964eb88SJed Brown                                        0,
5973964eb88SJed Brown                                        0,
5983964eb88SJed Brown                                        0,
5993964eb88SJed Brown                                /*124*/ 0,
6003964eb88SJed Brown                                        0,
6013964eb88SJed Brown                                        0,
6023964eb88SJed Brown                                        0,
6037dae84e0SHong Zhang                                        MatCreateSubMatricesMPI_MPIAdj,
6043964eb88SJed Brown                                /*129*/ 0,
6053964eb88SJed Brown                                        0,
6063964eb88SJed Brown                                        0,
6073964eb88SJed Brown                                        0,
6083964eb88SJed Brown                                        0,
6093964eb88SJed Brown                                /*134*/ 0,
6103964eb88SJed Brown                                        0,
6113964eb88SJed Brown                                        0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                        0,
6143964eb88SJed Brown                                /*139*/ 0,
615f9426fe0SMark Adams                                        0,
6163964eb88SJed Brown                                        0
6173964eb88SJed Brown };
618b97cf49bSBarry Smith 
619f7a08781SBarry Smith static PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
620a23d5eceSKris Buschelman {
621a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj*)B->data;
622*e895ccc0SFande Kong   PetscBool       useedgeweights;
623dfbe8321SBarry Smith   PetscErrorCode ierr;
6242515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
625b24ad042SBarry Smith   PetscInt ii;
626a23d5eceSKris Buschelman #endif
627a23d5eceSKris Buschelman 
628a23d5eceSKris Buschelman   PetscFunctionBegin;
62926283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
63026283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
631*e895ccc0SFande Kong   if (values) useedgeweights = PETSC_TRUE; else useedgeweights = PETSC_FALSE;
632*e895ccc0SFande Kong   /* Make everybody knows if they are using edge weights or not */
633*e895ccc0SFande Kong   ierr = MPIU_Allreduce((int*)&useedgeweights,(int*)&b->useedgeweights,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)B));CHKERRQ(ierr);
63458ec18a5SLisandro Dalcin 
6352515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
636e32f2f54SBarry 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]);
637d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
638e02043d6SBarry 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]);
639a23d5eceSKris Buschelman   }
640d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
641e02043d6SBarry 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]);
642a23d5eceSKris Buschelman   }
643a23d5eceSKris Buschelman #endif
64458ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
645a23d5eceSKris Buschelman 
646a23d5eceSKris Buschelman   b->j      = j;
647a23d5eceSKris Buschelman   b->i      = i;
648a23d5eceSKris Buschelman   b->values = values;
649a23d5eceSKris Buschelman 
650d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
651a23d5eceSKris Buschelman   b->diag      = 0;
652a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
653a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
654a23d5eceSKris Buschelman 
655a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
656a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
657a23d5eceSKris Buschelman   PetscFunctionReturn(0);
658a23d5eceSKris Buschelman }
659b97cf49bSBarry Smith 
660f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B)
6619a3665c6SJed Brown {
6629a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
6639a3665c6SJed Brown   PetscErrorCode ierr;
6649a3665c6SJed Brown   const PetscInt *ranges;
6659a3665c6SJed Brown   MPI_Comm       acomm,bcomm;
6669a3665c6SJed Brown   MPI_Group      agroup,bgroup;
6679a3665c6SJed Brown   PetscMPIInt    i,rank,size,nranks,*ranks;
6689a3665c6SJed Brown 
6699a3665c6SJed Brown   PetscFunctionBegin;
6700298fd71SBarry Smith   *B    = NULL;
671ce94432eSBarry Smith   ierr  = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr);
6729a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&size);CHKERRQ(ierr);
6739a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr);
6749a3665c6SJed Brown   ierr  = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
6759a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
6769a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) nranks++;
6779a3665c6SJed Brown   }
6789a3665c6SJed Brown   if (nranks == size) {         /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
6799a3665c6SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
6809a3665c6SJed Brown     *B   = A;
6819a3665c6SJed Brown     PetscFunctionReturn(0);
6829a3665c6SJed Brown   }
6839a3665c6SJed Brown 
684785e854fSJed Brown   ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr);
6859a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
6869a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i;
6879a3665c6SJed Brown   }
6889a3665c6SJed Brown   ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr);
6899a3665c6SJed Brown   ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr);
6909a3665c6SJed Brown   ierr = PetscFree(ranks);CHKERRQ(ierr);
6919a3665c6SJed Brown   ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr);
6929a3665c6SJed Brown   ierr = MPI_Group_free(&agroup);CHKERRQ(ierr);
6939a3665c6SJed Brown   ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr);
6949a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
6959a3665c6SJed Brown     PetscInt   m,N;
6969a3665c6SJed Brown     Mat_MPIAdj *b;
6970298fd71SBarry Smith     ierr       = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
6980298fd71SBarry Smith     ierr       = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
6999a3665c6SJed Brown     ierr       = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr);
7009a3665c6SJed Brown     b          = (Mat_MPIAdj*)(*B)->data;
7019a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
7029a3665c6SJed Brown     ierr       = MPI_Comm_free(&bcomm);CHKERRQ(ierr);
7039a3665c6SJed Brown   }
7049a3665c6SJed Brown   PetscFunctionReturn(0);
7059a3665c6SJed Brown }
7069a3665c6SJed Brown 
707b44e7856SBarry Smith PetscErrorCode  MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B)
708b44e7856SBarry Smith {
709b44e7856SBarry Smith   PetscErrorCode ierr;
710b44e7856SBarry Smith   PetscInt       M,N,*II,*J,NZ,nz,m,nzstart,i;
711b44e7856SBarry Smith   PetscInt       *Values = NULL;
712b44e7856SBarry Smith   Mat_MPIAdj     *adj = (Mat_MPIAdj*)A->data;
713b44e7856SBarry Smith   PetscMPIInt    mnz,mm,*allnz,*allm,size,*dispnz,*dispm;
714b44e7856SBarry Smith 
715b44e7856SBarry Smith   PetscFunctionBegin;
716b44e7856SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
717b44e7856SBarry Smith   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
718b44e7856SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
719b44e7856SBarry Smith   nz   = adj->nz;
720d4e69552SBarry Smith   if (adj->i[m] != nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nz %D not correct i[m] %d",nz,adj->i[m]);
721b44e7856SBarry Smith   ierr = MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
722d4e69552SBarry Smith 
723b44e7856SBarry Smith   ierr = PetscMPIIntCast(nz,&mnz);CHKERRQ(ierr);
724d4e69552SBarry Smith   ierr = PetscCalloc2(size,&allnz,size,&dispnz);CHKERRQ(ierr);
725d4e69552SBarry Smith   ierr = MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
726b44e7856SBarry Smith   dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1];
727b44e7856SBarry Smith   if (adj->values) {
728b44e7856SBarry Smith     ierr = PetscMalloc1(NZ,&Values);CHKERRQ(ierr);
729b44e7856SBarry Smith     ierr = MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
730b44e7856SBarry Smith   }
731b44e7856SBarry Smith   ierr = PetscMalloc1(NZ,&J);CHKERRQ(ierr);
732b44e7856SBarry Smith   ierr = MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
733b44e7856SBarry Smith   ierr = PetscFree2(allnz,dispnz);CHKERRQ(ierr);
734b44e7856SBarry Smith   ierr = MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
735b44e7856SBarry Smith   nzstart -= nz;
736b44e7856SBarry Smith   /* shift the i[] values so they will be correct after being received */
737b44e7856SBarry Smith   for (i=0; i<m; i++) adj->i[i] += nzstart;
738b44e7856SBarry Smith   ierr = PetscMalloc1(M+1,&II);CHKERRQ(ierr);
739b44e7856SBarry Smith   ierr = PetscMPIIntCast(m,&mm);CHKERRQ(ierr);
740b44e7856SBarry Smith   ierr = PetscMalloc2(size,&allm,size,&dispm);CHKERRQ(ierr);
741d4e69552SBarry Smith   ierr = MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
742b44e7856SBarry Smith   dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1];
743b44e7856SBarry Smith   ierr = MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
744b44e7856SBarry Smith   ierr = PetscFree2(allm,dispm);CHKERRQ(ierr);
745b44e7856SBarry Smith   II[M] = NZ;
746b44e7856SBarry Smith   /* shift the i[] values back */
747b44e7856SBarry Smith   for (i=0; i<m; i++) adj->i[i] -= nzstart;
748b44e7856SBarry Smith   ierr = MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B);CHKERRQ(ierr);
749b44e7856SBarry Smith   PetscFunctionReturn(0);
750b44e7856SBarry Smith }
751b44e7856SBarry Smith 
7529a3665c6SJed Brown /*@
7539a3665c6SJed Brown    MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows
7549a3665c6SJed Brown 
7559a3665c6SJed Brown    Collective
7569a3665c6SJed Brown 
7579a3665c6SJed Brown    Input Arguments:
7589a3665c6SJed Brown .  A - original MPIAdj matrix
7599a3665c6SJed Brown 
7609a3665c6SJed Brown    Output Arguments:
7610298fd71SBarry Smith .  B - matrix on subcommunicator, NULL on ranks that owned zero rows of A
7629a3665c6SJed Brown 
7639a3665c6SJed Brown    Level: developer
7649a3665c6SJed Brown 
7659a3665c6SJed Brown    Note:
7669a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
7679a3665c6SJed Brown 
7689a3665c6SJed Brown    The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed.
7699a3665c6SJed Brown 
7709a3665c6SJed Brown .seealso: MatCreateMPIAdj()
7719a3665c6SJed Brown @*/
7729a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B)
7739a3665c6SJed Brown {
7749a3665c6SJed Brown   PetscErrorCode ierr;
7759a3665c6SJed Brown 
7769a3665c6SJed Brown   PetscFunctionBegin;
7779a3665c6SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
7789a3665c6SJed Brown   ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
7799a3665c6SJed Brown   PetscFunctionReturn(0);
7809a3665c6SJed Brown }
7819a3665c6SJed Brown 
7820bad9183SKris Buschelman /*MC
783fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
7840bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
7850bad9183SKris Buschelman 
7860bad9183SKris Buschelman   Level: beginner
7870bad9183SKris Buschelman 
7880bad9183SKris Buschelman .seealso: MatCreateMPIAdj
7890bad9183SKris Buschelman M*/
7900bad9183SKris Buschelman 
7918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
792273d9f13SBarry Smith {
793273d9f13SBarry Smith   Mat_MPIAdj     *b;
7946849ba73SBarry Smith   PetscErrorCode ierr;
795273d9f13SBarry Smith 
796273d9f13SBarry Smith   PetscFunctionBegin;
797b00a9115SJed Brown   ierr         = PetscNewLog(B,&b);CHKERRQ(ierr);
798b0a32e0cSBarry Smith   B->data      = (void*)b;
799273d9f13SBarry Smith   ierr         = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
800273d9f13SBarry Smith   B->assembled = PETSC_FALSE;
801273d9f13SBarry Smith 
802bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
803bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr);
804b44e7856SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj);CHKERRQ(ierr);
80517667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
806273d9f13SBarry Smith   PetscFunctionReturn(0);
807273d9f13SBarry Smith }
808273d9f13SBarry Smith 
809a23d5eceSKris Buschelman /*@C
810b44e7856SBarry Smith    MatMPIAdjToSeq - Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential preconditioners)
811b44e7856SBarry Smith 
812b44e7856SBarry Smith    Logically Collective on MPI_Comm
813b44e7856SBarry Smith 
814b44e7856SBarry Smith    Input Parameter:
815b44e7856SBarry Smith .  A - the matrix
816b44e7856SBarry Smith 
817b44e7856SBarry Smith    Output Parameter:
818b44e7856SBarry Smith .  B - the same matrix on all processes
819b44e7856SBarry Smith 
820b44e7856SBarry Smith    Level: intermediate
821b44e7856SBarry Smith 
822b44e7856SBarry Smith .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
823b44e7856SBarry Smith @*/
824b44e7856SBarry Smith PetscErrorCode  MatMPIAdjToSeq(Mat A,Mat *B)
825b44e7856SBarry Smith {
826b44e7856SBarry Smith   PetscErrorCode ierr;
827b44e7856SBarry Smith 
828b44e7856SBarry Smith   PetscFunctionBegin;
829b44e7856SBarry Smith   ierr = PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
830b44e7856SBarry Smith   PetscFunctionReturn(0);
831b44e7856SBarry Smith }
832b44e7856SBarry Smith 
833b44e7856SBarry Smith /*@C
834a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
835a23d5eceSKris Buschelman 
8363f9fe445SBarry Smith    Logically Collective on MPI_Comm
837a23d5eceSKris Buschelman 
838a23d5eceSKris Buschelman    Input Parameters:
839a23d5eceSKris Buschelman +  A - the matrix
840a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
841a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
842a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
843a23d5eceSKris Buschelman -  values - [optional] edge weights
844a23d5eceSKris Buschelman 
845a23d5eceSKris Buschelman    Level: intermediate
846a23d5eceSKris Buschelman 
847a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
848a23d5eceSKris Buschelman @*/
8497087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
850273d9f13SBarry Smith {
8514ac538c5SBarry Smith   PetscErrorCode ierr;
852273d9f13SBarry Smith 
853273d9f13SBarry Smith   PetscFunctionBegin;
8544ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
855273d9f13SBarry Smith   PetscFunctionReturn(0);
856273d9f13SBarry Smith }
857273d9f13SBarry Smith 
858b97cf49bSBarry Smith /*@C
8593eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
860b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
861b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
862b97cf49bSBarry Smith 
863ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
864ef5ee4d1SLois Curfman McInnes 
865b97cf49bSBarry Smith    Input Parameters:
866c2e958feSBarry Smith +  comm - MPI communicator
8670752156aSBarry Smith .  m - number of local rows
86858ec18a5SLisandro Dalcin .  N - number of global columns
869b97cf49bSBarry Smith .  i - the indices into j for the start of each row
870329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
871ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
872329f5518SBarry Smith -  values -[optional] edge weights
873b97cf49bSBarry Smith 
874b97cf49bSBarry Smith    Output Parameter:
875b97cf49bSBarry Smith .  A - the matrix
876b97cf49bSBarry Smith 
87715091d37SBarry Smith    Level: intermediate
87815091d37SBarry Smith 
87995452b02SPatrick Sanan    Notes:
88095452b02SPatrick Sanan     This matrix object does not support most matrix operations, include
8814bc6d8c0SBarry Smith    MatSetValues().
882c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
883ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
8841198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
885327e1a73SBarry Smith    Should not include the matrix diagonals.
886b97cf49bSBarry Smith 
887e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
888ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
889ededeb1aSvictorle 
890ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
891b97cf49bSBarry Smith 
892e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
893b97cf49bSBarry Smith @*/
8947087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
895b97cf49bSBarry Smith {
896dfbe8321SBarry Smith   PetscErrorCode ierr;
897b97cf49bSBarry Smith 
898433994e6SBarry Smith   PetscFunctionBegin;
899f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
90058ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
901273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
902273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
9033a40ed3dSBarry Smith   PetscFunctionReturn(0);
904b97cf49bSBarry Smith }
905b97cf49bSBarry Smith 
906c2e958feSBarry Smith 
907433994e6SBarry Smith 
908433994e6SBarry Smith 
909433994e6SBarry Smith 
910433994e6SBarry Smith 
911433994e6SBarry Smith 
912