xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 0170919e3807ff82287a31d6cf5edc59527b75e2)
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);
18762795ce3SStefano Zampini #if defined(PETSC_USE_DEBUG)
18862795ce3SStefano Zampini         if (loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %D",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++) {
252*0170919eSFande Kong         if (a->values) {
253*0170919eSFande Kong           ierr = PetscViewerASCIISynchronizedPrintf(viewer," (%D, %D) ",a->j[j], a->values[j]);CHKERRQ(ierr);
254*0170919eSFande Kong         } else {
25577431f27SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr);
2560752156aSBarry Smith         }
257*0170919eSFande Kong       }
258b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
259b97cf49bSBarry Smith     }
260d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
261b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
262c5e4d11fSDmitry Karpeev     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
2637b23a99aSBarry Smith   }
2643a40ed3dSBarry Smith   PetscFunctionReturn(0);
265b97cf49bSBarry Smith }
266b97cf49bSBarry Smith 
26740244768SBarry Smith static PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer)
268b97cf49bSBarry Smith {
269dfbe8321SBarry Smith   PetscErrorCode ierr;
270ace3abfcSBarry Smith   PetscBool      iascii;
271b97cf49bSBarry Smith 
272433994e6SBarry Smith   PetscFunctionBegin;
273251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
27432077d6dSBarry Smith   if (iascii) {
2753eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
276b97cf49bSBarry Smith   }
2773a40ed3dSBarry Smith   PetscFunctionReturn(0);
278b97cf49bSBarry Smith }
279b97cf49bSBarry Smith 
28040244768SBarry Smith static PetscErrorCode MatDestroy_MPIAdj(Mat mat)
281b97cf49bSBarry Smith {
2823eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)mat->data;
283dfbe8321SBarry Smith   PetscErrorCode ierr;
28494d884c6SBarry Smith 
28594d884c6SBarry Smith   PetscFunctionBegin;
286aa482453SBarry Smith #if defined(PETSC_USE_LOG)
287d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz);
288b97cf49bSBarry Smith #endif
28905b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
2900400557aSBarry Smith   if (a->freeaij) {
29114196391SBarry Smith     if (a->freeaijwithfree) {
29214196391SBarry Smith       if (a->i) free(a->i);
29314196391SBarry Smith       if (a->j) free(a->j);
29414196391SBarry Smith     } else {
295606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
296606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
29705b42c5fSBarry Smith       ierr = PetscFree(a->values);CHKERRQ(ierr);
29814196391SBarry Smith     }
2990400557aSBarry Smith   }
3002b1d8763SJed Brown   ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
301bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
302dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
303bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr);
304bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr);
3053a40ed3dSBarry Smith   PetscFunctionReturn(0);
306b97cf49bSBarry Smith }
307b97cf49bSBarry Smith 
30840244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg)
309b97cf49bSBarry Smith {
3103eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
31163ba0a88SBarry Smith   PetscErrorCode ierr;
312b97cf49bSBarry Smith 
313433994e6SBarry Smith   PetscFunctionBegin;
31412c028f9SKris Buschelman   switch (op) {
31577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
31612c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
3179a4540c5SBarry Smith   case MAT_HERMITIAN:
3184e0d8c25SBarry Smith     a->symmetric = flg;
3199a4540c5SBarry Smith     break;
3209a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
3219a4540c5SBarry Smith     break;
32212c028f9SKris Buschelman   default:
323290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
32412c028f9SKris Buschelman     break;
325b97cf49bSBarry Smith   }
3263a40ed3dSBarry Smith   PetscFunctionReturn(0);
327b97cf49bSBarry Smith }
328b97cf49bSBarry Smith 
32940244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
330b97cf49bSBarry Smith {
3313eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
3322b1d8763SJed Brown   PetscErrorCode ierr;
333b97cf49bSBarry Smith 
334433994e6SBarry Smith   PetscFunctionBegin;
335d0f46423SBarry Smith   row -= A->rmap->rstart;
336e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
337b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
3382b1d8763SJed Brown   if (v) {
3392b1d8763SJed Brown     PetscInt j;
3402b1d8763SJed Brown     if (a->rowvalues_alloc < *nz) {
3412b1d8763SJed Brown       ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
3422b1d8763SJed Brown       a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz);
343785e854fSJed Brown       ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr);
344b97cf49bSBarry Smith     }
34591f8cafeSFande Kong     for (j=0; j<*nz; j++){
34691f8cafeSFande Kong       a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0;
34791f8cafeSFande Kong     }
3482b1d8763SJed Brown     *v = (*nz) ? a->rowvalues : NULL;
349b97cf49bSBarry Smith   }
35092b6f2a0SJed Brown   if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
3513a40ed3dSBarry Smith   PetscFunctionReturn(0);
352b97cf49bSBarry Smith }
353b97cf49bSBarry Smith 
35440244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
355b97cf49bSBarry Smith {
356433994e6SBarry Smith   PetscFunctionBegin;
3573a40ed3dSBarry Smith   PetscFunctionReturn(0);
358b97cf49bSBarry Smith }
359b97cf49bSBarry Smith 
36040244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
361b97cf49bSBarry Smith {
3623eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data;
363dfbe8321SBarry Smith   PetscErrorCode ierr;
364ace3abfcSBarry Smith   PetscBool      flag;
365b97cf49bSBarry Smith 
366433994e6SBarry Smith   PetscFunctionBegin;
367b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
368d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
3690f5bd95cSBarry Smith     flag = PETSC_FALSE;
370b97cf49bSBarry Smith   }
371b97cf49bSBarry Smith 
372b97cf49bSBarry Smith   /* if the a->i are the same */
373d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
374b97cf49bSBarry Smith 
375b97cf49bSBarry Smith   /* if a->j are the same */
376b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
377b97cf49bSBarry Smith 
378b2566f29SBarry Smith   ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
3793a40ed3dSBarry Smith   PetscFunctionReturn(0);
380b97cf49bSBarry Smith }
381b97cf49bSBarry Smith 
38240244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
3836cd91f64SBarry Smith {
384b24ad042SBarry Smith   PetscInt   i;
3856cd91f64SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
3861a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
3876cd91f64SBarry Smith 
3886cd91f64SBarry Smith   PetscFunctionBegin;
389d0f46423SBarry Smith   *m    = A->rmap->n;
3906cd91f64SBarry Smith   *ia   = a->i;
3916cd91f64SBarry Smith   *ja   = a->j;
3926cd91f64SBarry Smith   *done = PETSC_TRUE;
393d42735a1SBarry Smith   if (oshift) {
394d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
395d42735a1SBarry Smith       (*ja)[i]++;
396d42735a1SBarry Smith     }
397d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
398d42735a1SBarry Smith   }
399d42735a1SBarry Smith   PetscFunctionReturn(0);
400d42735a1SBarry Smith }
401d42735a1SBarry Smith 
40240244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
403d42735a1SBarry Smith {
404b24ad042SBarry Smith   PetscInt   i;
405d42735a1SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
4061a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
407d42735a1SBarry Smith 
408d42735a1SBarry Smith   PetscFunctionBegin;
40970d8d27cSBarry Smith   if (ia && a->i != *ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
41070d8d27cSBarry Smith   if (ja && a->j != *ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
411d42735a1SBarry Smith   if (oshift) {
41270d8d27cSBarry Smith     if (!ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument");
41370d8d27cSBarry Smith     if (!ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument");
414d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
415d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
416d42735a1SBarry Smith       (*ja)[i]--;
417d42735a1SBarry Smith     }
418d42735a1SBarry Smith   }
4196cd91f64SBarry Smith   PetscFunctionReturn(0);
4206cd91f64SBarry Smith }
421b97cf49bSBarry Smith 
42219fd82e9SBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
42317667f90SBarry Smith {
42417667f90SBarry Smith   Mat               B;
42517667f90SBarry Smith   PetscErrorCode    ierr;
42617667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
42717667f90SBarry Smith   const PetscInt    *rj;
42817667f90SBarry Smith   const PetscScalar *ra;
42917667f90SBarry Smith   MPI_Comm          comm;
43017667f90SBarry Smith 
43117667f90SBarry Smith   PetscFunctionBegin;
4320298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
4330298fd71SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
4340298fd71SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
43517667f90SBarry Smith 
43617667f90SBarry Smith   /* count the number of nonzeros per row */
43717667f90SBarry Smith   for (i=0; i<m; i++) {
4380298fd71SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
4395ee9ba1cSJed Brown     for (j=0; j<len; j++) {
4405ee9ba1cSJed Brown       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
4415ee9ba1cSJed Brown     }
44217667f90SBarry Smith     nzeros += len;
4430f2063bfSTobin Isaac     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
44417667f90SBarry Smith   }
44517667f90SBarry Smith 
44617667f90SBarry Smith   /* malloc space for nonzeros */
447854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr);
448854ce69bSBarry Smith   ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr);
449854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr);
45017667f90SBarry Smith 
45117667f90SBarry Smith   nzeros = 0;
45217667f90SBarry Smith   ia[0]  = 0;
45317667f90SBarry Smith   for (i=0; i<m; i++) {
45417667f90SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
45517667f90SBarry Smith     cnt  = 0;
45617667f90SBarry Smith     for (j=0; j<len; j++) {
4575ee9ba1cSJed Brown       if (rj[j] != i+rstart) { /* if not diagonal */
45817667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
45917667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
46017667f90SBarry Smith       }
4615ee9ba1cSJed Brown     }
46217667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
46317667f90SBarry Smith     nzeros += cnt;
46417667f90SBarry Smith     ia[i+1] = nzeros;
46517667f90SBarry Smith   }
46617667f90SBarry Smith 
46717667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
46817667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
46958ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
47017667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
47117667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
47217667f90SBarry Smith 
473511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
47428be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
47517667f90SBarry Smith   } else {
47617667f90SBarry Smith     *newmat = B;
47717667f90SBarry Smith   }
47817667f90SBarry Smith   PetscFunctionReturn(0);
47917667f90SBarry Smith }
48017667f90SBarry Smith 
481b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
4820b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
4833eda8832SBarry Smith                                        MatGetRow_MPIAdj,
4843eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
485b97cf49bSBarry Smith                                        0,
48697304618SKris Buschelman                                 /* 4*/ 0,
487b97cf49bSBarry Smith                                        0,
488b97cf49bSBarry Smith                                        0,
489b97cf49bSBarry Smith                                        0,
4900b3b1487SBarry Smith                                        0,
4910b3b1487SBarry Smith                                        0,
49297304618SKris Buschelman                                 /*10*/ 0,
4930b3b1487SBarry Smith                                        0,
4940b3b1487SBarry Smith                                        0,
4950b3b1487SBarry Smith                                        0,
4960b3b1487SBarry Smith                                        0,
49797304618SKris Buschelman                                 /*15*/ 0,
4983eda8832SBarry Smith                                        MatEqual_MPIAdj,
4990b3b1487SBarry Smith                                        0,
5000b3b1487SBarry Smith                                        0,
5010b3b1487SBarry Smith                                        0,
50297304618SKris Buschelman                                 /*20*/ 0,
5030b3b1487SBarry Smith                                        0,
5043eda8832SBarry Smith                                        MatSetOption_MPIAdj,
5050b3b1487SBarry Smith                                        0,
506d519adbfSMatthew Knepley                                 /*24*/ 0,
5070b3b1487SBarry Smith                                        0,
5080b3b1487SBarry Smith                                        0,
5090b3b1487SBarry Smith                                        0,
5100b3b1487SBarry Smith                                        0,
511d519adbfSMatthew Knepley                                 /*29*/ 0,
5120b3b1487SBarry Smith                                        0,
513273d9f13SBarry Smith                                        0,
514273d9f13SBarry Smith                                        0,
5150b3b1487SBarry Smith                                        0,
516d519adbfSMatthew Knepley                                 /*34*/ 0,
5170b3b1487SBarry Smith                                        0,
5180b3b1487SBarry Smith                                        0,
5190b3b1487SBarry Smith                                        0,
5200b3b1487SBarry Smith                                        0,
521d519adbfSMatthew Knepley                                 /*39*/ 0,
5227dae84e0SHong Zhang                                        MatCreateSubMatrices_MPIAdj,
5230b3b1487SBarry Smith                                        0,
5240b3b1487SBarry Smith                                        0,
5250b3b1487SBarry Smith                                        0,
526d519adbfSMatthew Knepley                                 /*44*/ 0,
5270b3b1487SBarry Smith                                        0,
5287d68702bSBarry Smith                                        MatShift_Basic,
5290b3b1487SBarry Smith                                        0,
5300b3b1487SBarry Smith                                        0,
531d519adbfSMatthew Knepley                                 /*49*/ 0,
5326cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
533d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
534b97cf49bSBarry Smith                                        0,
535b97cf49bSBarry Smith                                        0,
536d519adbfSMatthew Knepley                                 /*54*/ 0,
537b97cf49bSBarry Smith                                        0,
5380752156aSBarry Smith                                        0,
5390752156aSBarry Smith                                        0,
5400b3b1487SBarry Smith                                        0,
541d519adbfSMatthew Knepley                                 /*59*/ 0,
542b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
543b9b97703SBarry Smith                                        MatView_MPIAdj,
54417667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
54597304618SKris Buschelman                                        0,
546d519adbfSMatthew Knepley                                 /*64*/ 0,
54797304618SKris Buschelman                                        0,
54897304618SKris Buschelman                                        0,
54997304618SKris Buschelman                                        0,
55097304618SKris Buschelman                                        0,
551d519adbfSMatthew Knepley                                 /*69*/ 0,
55297304618SKris Buschelman                                        0,
55397304618SKris Buschelman                                        0,
55497304618SKris Buschelman                                        0,
55597304618SKris Buschelman                                        0,
556d519adbfSMatthew Knepley                                 /*74*/ 0,
55797304618SKris Buschelman                                        0,
55897304618SKris Buschelman                                        0,
55997304618SKris Buschelman                                        0,
56097304618SKris Buschelman                                        0,
561d519adbfSMatthew Knepley                                 /*79*/ 0,
56297304618SKris Buschelman                                        0,
56397304618SKris Buschelman                                        0,
56497304618SKris Buschelman                                        0,
565865e5f61SKris Buschelman                                        0,
566d519adbfSMatthew Knepley                                 /*84*/ 0,
567865e5f61SKris Buschelman                                        0,
568865e5f61SKris Buschelman                                        0,
569865e5f61SKris Buschelman                                        0,
570865e5f61SKris Buschelman                                        0,
571d519adbfSMatthew Knepley                                 /*89*/ 0,
572865e5f61SKris Buschelman                                        0,
573865e5f61SKris Buschelman                                        0,
574865e5f61SKris Buschelman                                        0,
575865e5f61SKris Buschelman                                        0,
576d519adbfSMatthew Knepley                                 /*94*/ 0,
577865e5f61SKris Buschelman                                        0,
578865e5f61SKris Buschelman                                        0,
5793964eb88SJed Brown                                        0,
5803964eb88SJed Brown                                        0,
5813964eb88SJed Brown                                 /*99*/ 0,
5823964eb88SJed Brown                                        0,
5833964eb88SJed Brown                                        0,
5843964eb88SJed Brown                                        0,
5853964eb88SJed Brown                                        0,
5863964eb88SJed Brown                                /*104*/ 0,
5873964eb88SJed Brown                                        0,
5883964eb88SJed Brown                                        0,
5893964eb88SJed Brown                                        0,
5903964eb88SJed Brown                                        0,
5913964eb88SJed Brown                                /*109*/ 0,
5923964eb88SJed Brown                                        0,
5933964eb88SJed Brown                                        0,
5943964eb88SJed Brown                                        0,
5953964eb88SJed Brown                                        0,
5963964eb88SJed Brown                                /*114*/ 0,
5973964eb88SJed Brown                                        0,
5983964eb88SJed Brown                                        0,
5993964eb88SJed Brown                                        0,
6003964eb88SJed Brown                                        0,
6013964eb88SJed Brown                                /*119*/ 0,
6023964eb88SJed Brown                                        0,
6033964eb88SJed Brown                                        0,
6043964eb88SJed Brown                                        0,
6053964eb88SJed Brown                                        0,
6063964eb88SJed Brown                                /*124*/ 0,
6073964eb88SJed Brown                                        0,
6083964eb88SJed Brown                                        0,
6093964eb88SJed Brown                                        0,
6107dae84e0SHong Zhang                                        MatCreateSubMatricesMPI_MPIAdj,
6113964eb88SJed Brown                                /*129*/ 0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                        0,
6143964eb88SJed Brown                                        0,
6153964eb88SJed Brown                                        0,
6163964eb88SJed Brown                                /*134*/ 0,
6173964eb88SJed Brown                                        0,
6183964eb88SJed Brown                                        0,
6193964eb88SJed Brown                                        0,
6203964eb88SJed Brown                                        0,
6213964eb88SJed Brown                                /*139*/ 0,
622f9426fe0SMark Adams                                        0,
6233964eb88SJed Brown                                        0
6243964eb88SJed Brown };
625b97cf49bSBarry Smith 
626f7a08781SBarry Smith static PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
627a23d5eceSKris Buschelman {
628a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj*)B->data;
629dfbe8321SBarry Smith   PetscErrorCode ierr;
6302515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
631b24ad042SBarry Smith   PetscInt ii;
632a23d5eceSKris Buschelman #endif
633a23d5eceSKris Buschelman 
634a23d5eceSKris Buschelman   PetscFunctionBegin;
63526283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
63626283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
63758ec18a5SLisandro Dalcin 
6382515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
639e32f2f54SBarry 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]);
640d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
641e02043d6SBarry 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]);
642a23d5eceSKris Buschelman   }
643d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
644e02043d6SBarry 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]);
645a23d5eceSKris Buschelman   }
646a23d5eceSKris Buschelman #endif
64758ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
648a23d5eceSKris Buschelman 
649a23d5eceSKris Buschelman   b->j      = j;
650a23d5eceSKris Buschelman   b->i      = i;
651a23d5eceSKris Buschelman   b->values = values;
652a23d5eceSKris Buschelman 
653d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
654a23d5eceSKris Buschelman   b->diag      = 0;
655a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
656a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
657a23d5eceSKris Buschelman 
658a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
659a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
660a23d5eceSKris Buschelman   PetscFunctionReturn(0);
661a23d5eceSKris Buschelman }
662b97cf49bSBarry Smith 
663f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B)
6649a3665c6SJed Brown {
6659a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
6669a3665c6SJed Brown   PetscErrorCode ierr;
6679a3665c6SJed Brown   const PetscInt *ranges;
6689a3665c6SJed Brown   MPI_Comm       acomm,bcomm;
6699a3665c6SJed Brown   MPI_Group      agroup,bgroup;
6709a3665c6SJed Brown   PetscMPIInt    i,rank,size,nranks,*ranks;
6719a3665c6SJed Brown 
6729a3665c6SJed Brown   PetscFunctionBegin;
6730298fd71SBarry Smith   *B    = NULL;
674ce94432eSBarry Smith   ierr  = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr);
6759a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&size);CHKERRQ(ierr);
6769a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr);
6779a3665c6SJed Brown   ierr  = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
6789a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
6799a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) nranks++;
6809a3665c6SJed Brown   }
6819a3665c6SJed Brown   if (nranks == size) {         /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
6829a3665c6SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
6839a3665c6SJed Brown     *B   = A;
6849a3665c6SJed Brown     PetscFunctionReturn(0);
6859a3665c6SJed Brown   }
6869a3665c6SJed Brown 
687785e854fSJed Brown   ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr);
6889a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
6899a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i;
6909a3665c6SJed Brown   }
6919a3665c6SJed Brown   ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr);
6929a3665c6SJed Brown   ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr);
6939a3665c6SJed Brown   ierr = PetscFree(ranks);CHKERRQ(ierr);
6949a3665c6SJed Brown   ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr);
6959a3665c6SJed Brown   ierr = MPI_Group_free(&agroup);CHKERRQ(ierr);
6969a3665c6SJed Brown   ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr);
6979a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
6989a3665c6SJed Brown     PetscInt   m,N;
6999a3665c6SJed Brown     Mat_MPIAdj *b;
7000298fd71SBarry Smith     ierr       = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
7010298fd71SBarry Smith     ierr       = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
7029a3665c6SJed Brown     ierr       = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr);
7039a3665c6SJed Brown     b          = (Mat_MPIAdj*)(*B)->data;
7049a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
7059a3665c6SJed Brown     ierr       = MPI_Comm_free(&bcomm);CHKERRQ(ierr);
7069a3665c6SJed Brown   }
7079a3665c6SJed Brown   PetscFunctionReturn(0);
7089a3665c6SJed Brown }
7099a3665c6SJed Brown 
710b44e7856SBarry Smith PetscErrorCode  MatMPIAdjToSeq_MPIAdj(Mat A,Mat *B)
711b44e7856SBarry Smith {
712b44e7856SBarry Smith   PetscErrorCode ierr;
713b44e7856SBarry Smith   PetscInt       M,N,*II,*J,NZ,nz,m,nzstart,i;
714b44e7856SBarry Smith   PetscInt       *Values = NULL;
715b44e7856SBarry Smith   Mat_MPIAdj     *adj = (Mat_MPIAdj*)A->data;
716b44e7856SBarry Smith   PetscMPIInt    mnz,mm,*allnz,*allm,size,*dispnz,*dispm;
717b44e7856SBarry Smith 
718b44e7856SBarry Smith   PetscFunctionBegin;
719b44e7856SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
720b44e7856SBarry Smith   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
721b44e7856SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
722b44e7856SBarry Smith   nz   = adj->nz;
723d4e69552SBarry Smith   if (adj->i[m] != nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nz %D not correct i[m] %d",nz,adj->i[m]);
724b44e7856SBarry Smith   ierr = MPI_Allreduce(&nz,&NZ,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
725d4e69552SBarry Smith 
726b44e7856SBarry Smith   ierr = PetscMPIIntCast(nz,&mnz);CHKERRQ(ierr);
727d4e69552SBarry Smith   ierr = PetscCalloc2(size,&allnz,size,&dispnz);CHKERRQ(ierr);
728d4e69552SBarry Smith   ierr = MPI_Allgather(&mnz,1,MPI_INT,allnz,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
729b44e7856SBarry Smith   dispnz[0] = 0; for (i=1; i<size; i++) dispnz[i] = dispnz[i-1]+ allnz[i-1];
730b44e7856SBarry Smith   if (adj->values) {
731b44e7856SBarry Smith     ierr = PetscMalloc1(NZ,&Values);CHKERRQ(ierr);
732b44e7856SBarry Smith     ierr = MPI_Allgatherv(adj->values,mnz,MPIU_INT,Values,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
733b44e7856SBarry Smith   }
734b44e7856SBarry Smith   ierr = PetscMalloc1(NZ,&J);CHKERRQ(ierr);
735b44e7856SBarry Smith   ierr = MPI_Allgatherv(adj->j,mnz,MPIU_INT,J,allnz,dispnz,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
736b44e7856SBarry Smith   ierr = PetscFree2(allnz,dispnz);CHKERRQ(ierr);
737b44e7856SBarry Smith   ierr = MPI_Scan(&nz,&nzstart,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
738b44e7856SBarry Smith   nzstart -= nz;
739b44e7856SBarry Smith   /* shift the i[] values so they will be correct after being received */
740b44e7856SBarry Smith   for (i=0; i<m; i++) adj->i[i] += nzstart;
741b44e7856SBarry Smith   ierr = PetscMalloc1(M+1,&II);CHKERRQ(ierr);
742b44e7856SBarry Smith   ierr = PetscMPIIntCast(m,&mm);CHKERRQ(ierr);
743b44e7856SBarry Smith   ierr = PetscMalloc2(size,&allm,size,&dispm);CHKERRQ(ierr);
744d4e69552SBarry Smith   ierr = MPI_Allgather(&mm,1,MPI_INT,allm,1,MPI_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
745b44e7856SBarry Smith   dispm[0] = 0; for (i=1; i<size; i++) dispm[i] = dispm[i-1]+ allm[i-1];
746b44e7856SBarry Smith   ierr = MPI_Allgatherv(adj->i,mm,MPIU_INT,II,allm,dispm,MPIU_INT,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
747b44e7856SBarry Smith   ierr = PetscFree2(allm,dispm);CHKERRQ(ierr);
748b44e7856SBarry Smith   II[M] = NZ;
749b44e7856SBarry Smith   /* shift the i[] values back */
750b44e7856SBarry Smith   for (i=0; i<m; i++) adj->i[i] -= nzstart;
751b44e7856SBarry Smith   ierr = MatCreateMPIAdj(PETSC_COMM_SELF,M,N,II,J,Values,B);CHKERRQ(ierr);
752b44e7856SBarry Smith   PetscFunctionReturn(0);
753b44e7856SBarry Smith }
754b44e7856SBarry Smith 
7559a3665c6SJed Brown /*@
7569a3665c6SJed Brown    MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows
7579a3665c6SJed Brown 
7589a3665c6SJed Brown    Collective
7599a3665c6SJed Brown 
7609a3665c6SJed Brown    Input Arguments:
7619a3665c6SJed Brown .  A - original MPIAdj matrix
7629a3665c6SJed Brown 
7639a3665c6SJed Brown    Output Arguments:
7640298fd71SBarry Smith .  B - matrix on subcommunicator, NULL on ranks that owned zero rows of A
7659a3665c6SJed Brown 
7669a3665c6SJed Brown    Level: developer
7679a3665c6SJed Brown 
7689a3665c6SJed Brown    Note:
7699a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
7709a3665c6SJed Brown 
7719a3665c6SJed Brown    The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed.
7729a3665c6SJed Brown 
7739a3665c6SJed Brown .seealso: MatCreateMPIAdj()
7749a3665c6SJed Brown @*/
7759a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B)
7769a3665c6SJed Brown {
7779a3665c6SJed Brown   PetscErrorCode ierr;
7789a3665c6SJed Brown 
7799a3665c6SJed Brown   PetscFunctionBegin;
7809a3665c6SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
7819a3665c6SJed Brown   ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
7829a3665c6SJed Brown   PetscFunctionReturn(0);
7839a3665c6SJed Brown }
7849a3665c6SJed Brown 
7850bad9183SKris Buschelman /*MC
786fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
7870bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
7880bad9183SKris Buschelman 
7890bad9183SKris Buschelman   Level: beginner
7900bad9183SKris Buschelman 
7910bad9183SKris Buschelman .seealso: MatCreateMPIAdj
7920bad9183SKris Buschelman M*/
7930bad9183SKris Buschelman 
7948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
795273d9f13SBarry Smith {
796273d9f13SBarry Smith   Mat_MPIAdj     *b;
7976849ba73SBarry Smith   PetscErrorCode ierr;
798273d9f13SBarry Smith 
799273d9f13SBarry Smith   PetscFunctionBegin;
800b00a9115SJed Brown   ierr         = PetscNewLog(B,&b);CHKERRQ(ierr);
801b0a32e0cSBarry Smith   B->data      = (void*)b;
802273d9f13SBarry Smith   ierr         = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
803273d9f13SBarry Smith   B->assembled = PETSC_FALSE;
804273d9f13SBarry Smith 
805bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
806bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr);
807b44e7856SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjToSeq_C",MatMPIAdjToSeq_MPIAdj);CHKERRQ(ierr);
80817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
809273d9f13SBarry Smith   PetscFunctionReturn(0);
810273d9f13SBarry Smith }
811273d9f13SBarry Smith 
812a23d5eceSKris Buschelman /*@C
813b44e7856SBarry Smith    MatMPIAdjToSeq - Converts an parallel MPIAdj matrix to complete MPIAdj on each process (needed by sequential preconditioners)
814b44e7856SBarry Smith 
815b44e7856SBarry Smith    Logically Collective on MPI_Comm
816b44e7856SBarry Smith 
817b44e7856SBarry Smith    Input Parameter:
818b44e7856SBarry Smith .  A - the matrix
819b44e7856SBarry Smith 
820b44e7856SBarry Smith    Output Parameter:
821b44e7856SBarry Smith .  B - the same matrix on all processes
822b44e7856SBarry Smith 
823b44e7856SBarry Smith    Level: intermediate
824b44e7856SBarry Smith 
825b44e7856SBarry Smith .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
826b44e7856SBarry Smith @*/
827b44e7856SBarry Smith PetscErrorCode  MatMPIAdjToSeq(Mat A,Mat *B)
828b44e7856SBarry Smith {
829b44e7856SBarry Smith   PetscErrorCode ierr;
830b44e7856SBarry Smith 
831b44e7856SBarry Smith   PetscFunctionBegin;
832b44e7856SBarry Smith   ierr = PetscUseMethod(A,"MatMPIAdjToSeq_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
833b44e7856SBarry Smith   PetscFunctionReturn(0);
834b44e7856SBarry Smith }
835b44e7856SBarry Smith 
836b44e7856SBarry Smith /*@C
837a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
838a23d5eceSKris Buschelman 
8393f9fe445SBarry Smith    Logically Collective on MPI_Comm
840a23d5eceSKris Buschelman 
841a23d5eceSKris Buschelman    Input Parameters:
842a23d5eceSKris Buschelman +  A - the matrix
843a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
844a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
845a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
846a23d5eceSKris Buschelman -  values - [optional] edge weights
847a23d5eceSKris Buschelman 
848a23d5eceSKris Buschelman    Level: intermediate
849a23d5eceSKris Buschelman 
850a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
851a23d5eceSKris Buschelman @*/
8527087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
853273d9f13SBarry Smith {
8544ac538c5SBarry Smith   PetscErrorCode ierr;
855273d9f13SBarry Smith 
856273d9f13SBarry Smith   PetscFunctionBegin;
8574ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
858273d9f13SBarry Smith   PetscFunctionReturn(0);
859273d9f13SBarry Smith }
860273d9f13SBarry Smith 
861b97cf49bSBarry Smith /*@C
8623eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
863b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
864b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
865b97cf49bSBarry Smith 
866ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
867ef5ee4d1SLois Curfman McInnes 
868b97cf49bSBarry Smith    Input Parameters:
869c2e958feSBarry Smith +  comm - MPI communicator
8700752156aSBarry Smith .  m - number of local rows
87158ec18a5SLisandro Dalcin .  N - number of global columns
872b97cf49bSBarry Smith .  i - the indices into j for the start of each row
873329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
874ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
875329f5518SBarry Smith -  values -[optional] edge weights
876b97cf49bSBarry Smith 
877b97cf49bSBarry Smith    Output Parameter:
878b97cf49bSBarry Smith .  A - the matrix
879b97cf49bSBarry Smith 
88015091d37SBarry Smith    Level: intermediate
88115091d37SBarry Smith 
88295452b02SPatrick Sanan    Notes:
88395452b02SPatrick Sanan     This matrix object does not support most matrix operations, include
8844bc6d8c0SBarry Smith    MatSetValues().
885c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
886ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
8871198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
888327e1a73SBarry Smith    Should not include the matrix diagonals.
889b97cf49bSBarry Smith 
890e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
891ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
892ededeb1aSvictorle 
893ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
894b97cf49bSBarry Smith 
895e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
896b97cf49bSBarry Smith @*/
8977087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
898b97cf49bSBarry Smith {
899dfbe8321SBarry Smith   PetscErrorCode ierr;
900b97cf49bSBarry Smith 
901433994e6SBarry Smith   PetscFunctionBegin;
902f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
90358ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
904273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
905273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
9063a40ed3dSBarry Smith   PetscFunctionReturn(0);
907b97cf49bSBarry Smith }
908b97cf49bSBarry Smith 
909c2e958feSBarry Smith 
910433994e6SBarry Smith 
911433994e6SBarry Smith 
912433994e6SBarry Smith 
913433994e6SBarry Smith 
914433994e6SBarry Smith 
915