xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 40244768c647318636de500975f889962f70aaa7)
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 
84a2ae208SSatish Balay #undef __FUNCT__
9*40244768SBarry Smith #define __FUNCT__ "MatGetSubMatrix_MPIAdj_data"
10*40244768SBarry Smith /*
11*40244768SBarry Smith  * The interface should be easy to use for both MatGetSubMatrix (parallel sub-matrix) and MatGetSubMatrices (sequential sub-matrices)
12*40244768SBarry Smith  * */
13*40244768SBarry Smith static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat adj,IS irows, IS icols, PetscInt **sadj_xadj,PetscInt **sadj_adjncy,PetscInt **sadj_values)
14*40244768SBarry Smith {
15*40244768SBarry Smith   PetscInt           nlrows_is,icols_n,i,j,nroots,nleaves,owner,rlocalindex,*ncols_send,*ncols_recv;
16*40244768SBarry Smith   PetscInt           nlrows_mat,*adjncy_recv,Ncols_recv,Ncols_send,*xadj_recv,*values_recv;
17*40244768SBarry Smith   PetscInt          *ncols_recv_offsets,loc,rnclos,*sadjncy,*sxadj,*svalues,isvalue;
18*40244768SBarry Smith   const PetscInt    *irows_indices,*icols_indices,*xadj, *adjncy;
19*40244768SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)adj->data;
20*40244768SBarry Smith   PetscLayout        rmap;
21*40244768SBarry Smith   MPI_Comm           comm;
22*40244768SBarry Smith   PetscSF            sf;
23*40244768SBarry Smith   PetscSFNode       *iremote;
24*40244768SBarry Smith   PetscBool          done;
25*40244768SBarry Smith   PetscErrorCode     ierr;
26*40244768SBarry Smith 
27*40244768SBarry Smith   PetscFunctionBegin;
28*40244768SBarry Smith   ierr = PetscObjectGetComm((PetscObject)adj,&comm);CHKERRQ(ierr);
29*40244768SBarry Smith   ierr = MatGetLayouts(adj,&rmap,PETSC_NULL);CHKERRQ(ierr);
30*40244768SBarry Smith   ierr = ISGetLocalSize(irows,&nlrows_is);CHKERRQ(ierr);
31*40244768SBarry Smith   ierr = ISGetIndices(irows,&irows_indices);CHKERRQ(ierr);
32*40244768SBarry Smith   ierr = PetscCalloc1(nlrows_is,&iremote);CHKERRQ(ierr);
33*40244768SBarry Smith   /* construct sf graph*/
34*40244768SBarry Smith   nleaves = nlrows_is;
35*40244768SBarry Smith   for (i=0; i<nlrows_is; i++){
36*40244768SBarry Smith     owner = -1;
37*40244768SBarry Smith     rlocalindex = -1;
38*40244768SBarry Smith     ierr = PetscLayoutFindOwnerIndex(rmap,irows_indices[i],&owner,&rlocalindex);CHKERRQ(ierr);
39*40244768SBarry Smith     iremote[i].rank  = owner;
40*40244768SBarry Smith     iremote[i].index = rlocalindex;
41*40244768SBarry Smith   }
42*40244768SBarry Smith   ierr = MatGetRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr);
43*40244768SBarry Smith   ierr = PetscCalloc4(nlrows_mat,&ncols_send,nlrows_is,&xadj_recv,nlrows_is+1,&ncols_recv_offsets,nlrows_is,&ncols_recv);CHKERRQ(ierr);
44*40244768SBarry Smith   nroots = nlrows_mat;
45*40244768SBarry Smith   for (i=0; i<nlrows_mat; i++){
46*40244768SBarry Smith     ncols_send[i] = xadj[i+1]-xadj[i];
47*40244768SBarry Smith   }
48*40244768SBarry Smith   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
49*40244768SBarry Smith   ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
50*40244768SBarry Smith   ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);
51*40244768SBarry Smith   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
52*40244768SBarry Smith   ierr = PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr);
53*40244768SBarry Smith   ierr = PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr);
54*40244768SBarry Smith   ierr = PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr);
55*40244768SBarry Smith   ierr = PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr);
56*40244768SBarry Smith   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
57*40244768SBarry Smith   Ncols_recv =0;
58*40244768SBarry Smith   for (i=0; i<nlrows_is; i++){
59*40244768SBarry Smith     Ncols_recv             += ncols_recv[i];
60*40244768SBarry Smith     ncols_recv_offsets[i+1] = ncols_recv[i]+ncols_recv_offsets[i];
61*40244768SBarry Smith   }
62*40244768SBarry Smith   Ncols_send = 0;
63*40244768SBarry Smith   for(i=0; i<nlrows_mat; i++){
64*40244768SBarry Smith     Ncols_send += ncols_send[i];
65*40244768SBarry Smith   }
66*40244768SBarry Smith   ierr = PetscCalloc1(Ncols_recv,&iremote);CHKERRQ(ierr);
67*40244768SBarry Smith   ierr = PetscCalloc1(Ncols_recv,&adjncy_recv);CHKERRQ(ierr);
68*40244768SBarry Smith   nleaves = Ncols_recv;
69*40244768SBarry Smith   Ncols_recv = 0;
70*40244768SBarry Smith   for (i=0; i<nlrows_is; i++){
71*40244768SBarry Smith     ierr = PetscLayoutFindOwner(rmap,irows_indices[i],&owner);CHKERRQ(ierr);
72*40244768SBarry Smith     for (j=0; j<ncols_recv[i]; j++){
73*40244768SBarry Smith       iremote[Ncols_recv].rank    = owner;
74*40244768SBarry Smith       iremote[Ncols_recv++].index = xadj_recv[i]+j;
75*40244768SBarry Smith     }
76*40244768SBarry Smith   }
77*40244768SBarry Smith   ierr = ISRestoreIndices(irows,&irows_indices);CHKERRQ(ierr);
78*40244768SBarry Smith   /*if we need to deal with edge weights ???*/
79*40244768SBarry Smith   if (a->values){isvalue=1;} else {isvalue=0;}
80*40244768SBarry Smith   /*involve a global communication */
81*40244768SBarry Smith   /*ierr = MPIU_Allreduce(&isvalue,&isvalue,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);*/
82*40244768SBarry Smith   if (isvalue){ierr = PetscCalloc1(Ncols_recv,&values_recv);CHKERRQ(ierr);}
83*40244768SBarry Smith   nroots = Ncols_send;
84*40244768SBarry Smith   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
85*40244768SBarry Smith   ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
86*40244768SBarry Smith   ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);
87*40244768SBarry Smith   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
88*40244768SBarry Smith   ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr);
89*40244768SBarry Smith   ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr);
90*40244768SBarry Smith   if (isvalue){
91*40244768SBarry Smith     ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr);
92*40244768SBarry Smith     ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr);
93*40244768SBarry Smith   }
94*40244768SBarry Smith   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
95*40244768SBarry Smith   ierr = MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr);
96*40244768SBarry Smith   ierr = ISGetLocalSize(icols,&icols_n);CHKERRQ(ierr);
97*40244768SBarry Smith   ierr = ISGetIndices(icols,&icols_indices);CHKERRQ(ierr);
98*40244768SBarry Smith   rnclos = 0;
99*40244768SBarry Smith   for (i=0; i<nlrows_is; i++){
100*40244768SBarry Smith     for (j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){
101*40244768SBarry Smith       ierr = PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc);CHKERRQ(ierr);
102*40244768SBarry Smith       if (loc<0){
103*40244768SBarry Smith         adjncy_recv[j] = -1;
104*40244768SBarry Smith         if (isvalue) values_recv[j] = -1;
105*40244768SBarry Smith         ncols_recv[i]--;
106*40244768SBarry Smith       } else {
107*40244768SBarry Smith         rnclos++;
108*40244768SBarry Smith       }
109*40244768SBarry Smith     }
110*40244768SBarry Smith   }
111*40244768SBarry Smith   ierr = ISRestoreIndices(icols,&icols_indices);CHKERRQ(ierr);
112*40244768SBarry Smith   ierr = PetscCalloc1(rnclos,&sadjncy);CHKERRQ(ierr);
113*40244768SBarry Smith   if (isvalue) {ierr = PetscCalloc1(rnclos,&svalues);CHKERRQ(ierr);}
114*40244768SBarry Smith   ierr = PetscCalloc1(nlrows_is+1,&sxadj);CHKERRQ(ierr);
115*40244768SBarry Smith   rnclos = 0;
116*40244768SBarry Smith   for(i=0; i<nlrows_is; i++){
117*40244768SBarry Smith     for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){
118*40244768SBarry Smith       if (adjncy_recv[j]<0) continue;
119*40244768SBarry Smith       sadjncy[rnclos] = adjncy_recv[j];
120*40244768SBarry Smith       if (isvalue) svalues[rnclos] = values_recv[j];
121*40244768SBarry Smith       rnclos++;
122*40244768SBarry Smith     }
123*40244768SBarry Smith   }
124*40244768SBarry Smith   for(i=0; i<nlrows_is; i++){
125*40244768SBarry Smith     sxadj[i+1] = sxadj[i]+ncols_recv[i];
126*40244768SBarry Smith   }
127*40244768SBarry Smith   if (sadj_xadj)  { *sadj_xadj = sxadj;} else    { ierr = PetscFree(sxadj);CHKERRQ(ierr);}
128*40244768SBarry Smith   if (sadj_adjncy){ *sadj_adjncy = sadjncy;} else { ierr = PetscFree(sadjncy);CHKERRQ(ierr);}
129*40244768SBarry Smith   if (sadj_values){
130*40244768SBarry Smith     if (isvalue) *sadj_values = svalues; else *sadj_values=0;
131*40244768SBarry Smith   } else {
132*40244768SBarry Smith     if (isvalue) {ierr = PetscFree(svalues);CHKERRQ(ierr);}
133*40244768SBarry Smith   }
134*40244768SBarry Smith   ierr = PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv);CHKERRQ(ierr);
135*40244768SBarry Smith   ierr = PetscFree(adjncy_recv);CHKERRQ(ierr);
136*40244768SBarry Smith   if (isvalue) {ierr = PetscFree(values_recv);CHKERRQ(ierr);}
137*40244768SBarry Smith   PetscFunctionReturn(0);
138*40244768SBarry Smith }
139*40244768SBarry Smith 
140*40244768SBarry Smith #undef __FUNCT__
141*40244768SBarry Smith #define __FUNCT__ "MatGetSubMatrices_MPIAdj_Private"
142*40244768SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIAdj_Private(Mat mat,PetscInt n,const IS irow[],const IS icol[],PetscBool subcomm,MatReuse scall,Mat *submat[])
143*40244768SBarry Smith {
144*40244768SBarry Smith   PetscInt           i,irow_n,icol_n,*sxadj,*sadjncy,*svalues;
145*40244768SBarry Smith   PetscInt          *indices,nindx,j,k,loc;
146*40244768SBarry Smith   PetscMPIInt        issame;
147*40244768SBarry Smith   const PetscInt    *irow_indices,*icol_indices;
148*40244768SBarry Smith   MPI_Comm           scomm_row,scomm_col,scomm_mat;
149*40244768SBarry Smith   PetscErrorCode     ierr;
150*40244768SBarry Smith 
151*40244768SBarry Smith   PetscFunctionBegin;
152*40244768SBarry Smith   nindx = 0;
153*40244768SBarry Smith   /*
154*40244768SBarry Smith    * Estimate a maximum number for allocating memory
155*40244768SBarry Smith    */
156*40244768SBarry Smith   for(i=0; i<n; i++){
157*40244768SBarry Smith     ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr);
158*40244768SBarry Smith     ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr);
159*40244768SBarry Smith     nindx = nindx>(irow_n+icol_n)? nindx:(irow_n+icol_n);
160*40244768SBarry Smith   }
161*40244768SBarry Smith   ierr = PetscCalloc1(nindx,&indices);CHKERRQ(ierr);
162*40244768SBarry Smith   /* construct a submat */
163*40244768SBarry Smith   for(i=0; i<n; i++){
164*40244768SBarry Smith     if (subcomm){
165*40244768SBarry Smith       ierr = PetscObjectGetComm((PetscObject)irow[i],&scomm_row);CHKERRQ(ierr);
166*40244768SBarry Smith       ierr = PetscObjectGetComm((PetscObject)icol[i],&scomm_col);CHKERRQ(ierr);
167*40244768SBarry Smith       ierr = MPI_Comm_compare(scomm_row,scomm_col,&issame);CHKERRQ(ierr);
168*40244768SBarry 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");
169*40244768SBarry Smith       ierr = MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame);CHKERRQ(ierr);
170*40244768SBarry 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");
171*40244768SBarry Smith     } else {
172*40244768SBarry Smith       scomm_row = PETSC_COMM_SELF;
173*40244768SBarry Smith     }
174*40244768SBarry Smith     /*get sub-matrix data*/
175*40244768SBarry Smith     sxadj=0; sadjncy=0; svalues=0;
176*40244768SBarry Smith     ierr = MatGetSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues);CHKERRQ(ierr);
177*40244768SBarry Smith     ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr);
178*40244768SBarry Smith     ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr);
179*40244768SBarry Smith     ierr = ISGetIndices(irow[i],&irow_indices);CHKERRQ(ierr);
180*40244768SBarry Smith     ierr = PetscMemcpy(indices,irow_indices,sizeof(PetscInt)*irow_n);CHKERRQ(ierr);
181*40244768SBarry Smith     ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr);
182*40244768SBarry Smith     ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr);
183*40244768SBarry Smith     ierr = PetscMemcpy(indices+irow_n,icol_indices,sizeof(PetscInt)*icol_n);CHKERRQ(ierr);
184*40244768SBarry Smith     ierr = ISRestoreIndices(icol[i],&icol_indices);CHKERRQ(ierr);
185*40244768SBarry Smith     nindx = irow_n+icol_n;
186*40244768SBarry Smith     ierr = PetscSortRemoveDupsInt(&nindx,indices);CHKERRQ(ierr);
187*40244768SBarry Smith     /* renumber columns */
188*40244768SBarry Smith     for(j=0; j<irow_n; j++){
189*40244768SBarry Smith       for(k=sxadj[j]; k<sxadj[j+1]; k++){
190*40244768SBarry Smith         ierr = PetscFindInt(sadjncy[k],nindx,indices,&loc);CHKERRQ(ierr);
191*40244768SBarry Smith #if PETSC_USE_DEBUG
192*40244768SBarry Smith         if (loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %d \n",sadjncy[k]);
193*40244768SBarry Smith #endif
194*40244768SBarry Smith         sadjncy[k] = loc;
195*40244768SBarry Smith       }
196*40244768SBarry Smith     }
197*40244768SBarry Smith     if (scall==MAT_INITIAL_MATRIX){
198*40244768SBarry Smith       ierr = MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr);
199*40244768SBarry Smith     } else {
200*40244768SBarry Smith        Mat                sadj = *(submat[i]);
201*40244768SBarry Smith        Mat_MPIAdj         *sa  = (Mat_MPIAdj*)((sadj)->data);
202*40244768SBarry Smith        ierr = PetscObjectGetComm((PetscObject)sadj,&scomm_mat);CHKERRQ(ierr);
203*40244768SBarry Smith        ierr = MPI_Comm_compare(scomm_row,scomm_mat,&issame);CHKERRQ(ierr);
204*40244768SBarry 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");
205*40244768SBarry Smith        ierr = PetscMemcpy(sa->i,sxadj,sizeof(PetscInt)*(irow_n+1));CHKERRQ(ierr);
206*40244768SBarry Smith        ierr = PetscMemcpy(sa->j,sadjncy,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);
207*40244768SBarry Smith        if (svalues){ierr = PetscMemcpy(sa->values,svalues,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);}
208*40244768SBarry Smith        ierr = PetscFree(sxadj);CHKERRQ(ierr);
209*40244768SBarry Smith        ierr = PetscFree(sadjncy);CHKERRQ(ierr);
210*40244768SBarry Smith        if (svalues) {ierr = PetscFree(svalues);CHKERRQ(ierr);}
211*40244768SBarry Smith     }
212*40244768SBarry Smith   }
213*40244768SBarry Smith   ierr = PetscFree(indices);CHKERRQ(ierr);
214*40244768SBarry Smith   PetscFunctionReturn(0);
215*40244768SBarry Smith }
216*40244768SBarry Smith 
217*40244768SBarry Smith #undef __FUNCT__
218*40244768SBarry Smith #define __FUNCT__ "MatGetSubMatricesMPI_MPIAdj"
219*40244768SBarry Smith static PetscErrorCode MatGetSubMatricesMPI_MPIAdj(Mat mat,PetscInt n, const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
220*40244768SBarry Smith {
221*40244768SBarry Smith   PetscErrorCode     ierr;
222*40244768SBarry Smith   /*get sub-matrices across a sub communicator */
223*40244768SBarry Smith   PetscFunctionBegin;
224*40244768SBarry Smith   ierr = MatGetSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat);CHKERRQ(ierr);
225*40244768SBarry Smith   PetscFunctionReturn(0);
226*40244768SBarry Smith }
227*40244768SBarry Smith 
228*40244768SBarry Smith #undef __FUNCT__
229*40244768SBarry Smith #define __FUNCT__ "MatGetSubMatrices_MPIAdj"
230*40244768SBarry Smith static PetscErrorCode MatGetSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
231*40244768SBarry Smith {
232*40244768SBarry Smith   PetscErrorCode     ierr;
233*40244768SBarry Smith 
234*40244768SBarry Smith   PetscFunctionBegin;
235*40244768SBarry Smith   /*get sub-matrices based on PETSC_COMM_SELF */
236*40244768SBarry Smith   ierr = MatGetSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat);CHKERRQ(ierr);
237*40244768SBarry Smith   PetscFunctionReturn(0);
238*40244768SBarry Smith }
239*40244768SBarry Smith 
240*40244768SBarry Smith #undef __FUNCT__
2414a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII"
242*40244768SBarry Smith static PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
243b97cf49bSBarry Smith {
2443eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
245dfbe8321SBarry Smith   PetscErrorCode    ierr;
246d0f46423SBarry Smith   PetscInt          i,j,m = A->rmap->n;
2472dcb1b2aSMatthew Knepley   const char        *name;
248ce1411ecSBarry Smith   PetscViewerFormat format;
249b97cf49bSBarry Smith 
250433994e6SBarry Smith   PetscFunctionBegin;
2513a7fca6bSBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
252b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
253fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
2543a40ed3dSBarry Smith     PetscFunctionReturn(0);
2556c4ed002SBarry Smith   } else if (format == PETSC_VIEWER_ASCII_MATLAB) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported");
2566c4ed002SBarry Smith   else {
257d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
258c5e4d11fSDmitry Karpeev     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
259b97cf49bSBarry Smith     for (i=0; i<m; i++) {
260d0f46423SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr);
261b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
26277431f27SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr);
2630752156aSBarry Smith       }
264b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
265b97cf49bSBarry Smith     }
266d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
267b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
268c5e4d11fSDmitry Karpeev     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
2697b23a99aSBarry Smith   }
2703a40ed3dSBarry Smith   PetscFunctionReturn(0);
271b97cf49bSBarry Smith }
272b97cf49bSBarry Smith 
2734a2ae208SSatish Balay #undef __FUNCT__
2744a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj"
275*40244768SBarry Smith static PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer)
276b97cf49bSBarry Smith {
277dfbe8321SBarry Smith   PetscErrorCode ierr;
278ace3abfcSBarry Smith   PetscBool      iascii;
279b97cf49bSBarry Smith 
280433994e6SBarry Smith   PetscFunctionBegin;
281251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
28232077d6dSBarry Smith   if (iascii) {
2833eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
284b97cf49bSBarry Smith   }
2853a40ed3dSBarry Smith   PetscFunctionReturn(0);
286b97cf49bSBarry Smith }
287b97cf49bSBarry Smith 
2884a2ae208SSatish Balay #undef __FUNCT__
2894a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj"
290*40244768SBarry Smith static PetscErrorCode MatDestroy_MPIAdj(Mat mat)
291b97cf49bSBarry Smith {
2923eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)mat->data;
293dfbe8321SBarry Smith   PetscErrorCode ierr;
29494d884c6SBarry Smith 
29594d884c6SBarry Smith   PetscFunctionBegin;
296aa482453SBarry Smith #if defined(PETSC_USE_LOG)
297d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz);
298b97cf49bSBarry Smith #endif
29905b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
3000400557aSBarry Smith   if (a->freeaij) {
30114196391SBarry Smith     if (a->freeaijwithfree) {
30214196391SBarry Smith       if (a->i) free(a->i);
30314196391SBarry Smith       if (a->j) free(a->j);
30414196391SBarry Smith     } else {
305606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
306606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
30705b42c5fSBarry Smith       ierr = PetscFree(a->values);CHKERRQ(ierr);
30814196391SBarry Smith     }
3090400557aSBarry Smith   }
3102b1d8763SJed Brown   ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
311bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
312dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
313bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr);
314bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr);
3153a40ed3dSBarry Smith   PetscFunctionReturn(0);
316b97cf49bSBarry Smith }
317b97cf49bSBarry Smith 
3184a2ae208SSatish Balay #undef __FUNCT__
3194a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
320*40244768SBarry Smith static PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg)
321b97cf49bSBarry Smith {
3223eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
32363ba0a88SBarry Smith   PetscErrorCode ierr;
324b97cf49bSBarry Smith 
325433994e6SBarry Smith   PetscFunctionBegin;
32612c028f9SKris Buschelman   switch (op) {
32777e54ba9SKris Buschelman   case MAT_SYMMETRIC:
32812c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
3299a4540c5SBarry Smith   case MAT_HERMITIAN:
3304e0d8c25SBarry Smith     a->symmetric = flg;
3319a4540c5SBarry Smith     break;
3329a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
3339a4540c5SBarry Smith     break;
33412c028f9SKris Buschelman   default:
335290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
33612c028f9SKris Buschelman     break;
337b97cf49bSBarry Smith   }
3383a40ed3dSBarry Smith   PetscFunctionReturn(0);
339b97cf49bSBarry Smith }
340b97cf49bSBarry Smith 
3414a2ae208SSatish Balay #undef __FUNCT__
3424a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
343*40244768SBarry Smith static PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
344b97cf49bSBarry Smith {
3453eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
3462b1d8763SJed Brown   PetscErrorCode ierr;
347b97cf49bSBarry Smith 
348433994e6SBarry Smith   PetscFunctionBegin;
349d0f46423SBarry Smith   row -= A->rmap->rstart;
350e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
351b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
3522b1d8763SJed Brown   if (v) {
3532b1d8763SJed Brown     PetscInt j;
3542b1d8763SJed Brown     if (a->rowvalues_alloc < *nz) {
3552b1d8763SJed Brown       ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
3562b1d8763SJed Brown       a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz);
357785e854fSJed Brown       ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr);
358b97cf49bSBarry Smith     }
35991f8cafeSFande Kong     for (j=0; j<*nz; j++){
36091f8cafeSFande Kong       a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0;
36191f8cafeSFande Kong     }
3622b1d8763SJed Brown     *v = (*nz) ? a->rowvalues : NULL;
363b97cf49bSBarry Smith   }
36492b6f2a0SJed Brown   if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
3653a40ed3dSBarry Smith   PetscFunctionReturn(0);
366b97cf49bSBarry Smith }
367b97cf49bSBarry Smith 
3684a2ae208SSatish Balay #undef __FUNCT__
3694a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
370*40244768SBarry Smith static PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
371b97cf49bSBarry Smith {
372433994e6SBarry Smith   PetscFunctionBegin;
3733a40ed3dSBarry Smith   PetscFunctionReturn(0);
374b97cf49bSBarry Smith }
375b97cf49bSBarry Smith 
3764a2ae208SSatish Balay #undef __FUNCT__
3774a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
378*40244768SBarry Smith static PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
379b97cf49bSBarry Smith {
3803eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data;
381dfbe8321SBarry Smith   PetscErrorCode ierr;
382ace3abfcSBarry Smith   PetscBool      flag;
383b97cf49bSBarry Smith 
384433994e6SBarry Smith   PetscFunctionBegin;
385b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
386d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
3870f5bd95cSBarry Smith     flag = PETSC_FALSE;
388b97cf49bSBarry Smith   }
389b97cf49bSBarry Smith 
390b97cf49bSBarry Smith   /* if the a->i are the same */
391d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
392b97cf49bSBarry Smith 
393b97cf49bSBarry Smith   /* if a->j are the same */
394b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
395b97cf49bSBarry Smith 
396b2566f29SBarry Smith   ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
3973a40ed3dSBarry Smith   PetscFunctionReturn(0);
398b97cf49bSBarry Smith }
399b97cf49bSBarry Smith 
4004a2ae208SSatish Balay #undef __FUNCT__
4014a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
402*40244768SBarry Smith static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
4036cd91f64SBarry Smith {
404b24ad042SBarry Smith   PetscInt   i;
4056cd91f64SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
4061a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
4076cd91f64SBarry Smith 
4086cd91f64SBarry Smith   PetscFunctionBegin;
409d0f46423SBarry Smith   *m    = A->rmap->n;
4106cd91f64SBarry Smith   *ia   = a->i;
4116cd91f64SBarry Smith   *ja   = a->j;
4126cd91f64SBarry Smith   *done = PETSC_TRUE;
413d42735a1SBarry Smith   if (oshift) {
414d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
415d42735a1SBarry Smith       (*ja)[i]++;
416d42735a1SBarry Smith     }
417d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
418d42735a1SBarry Smith   }
419d42735a1SBarry Smith   PetscFunctionReturn(0);
420d42735a1SBarry Smith }
421d42735a1SBarry Smith 
4224a2ae208SSatish Balay #undef __FUNCT__
4234a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
424*40244768SBarry Smith static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
425d42735a1SBarry Smith {
426b24ad042SBarry Smith   PetscInt   i;
427d42735a1SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
4281a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
429d42735a1SBarry Smith 
430d42735a1SBarry Smith   PetscFunctionBegin;
43170d8d27cSBarry Smith   if (ia && a->i != *ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
43270d8d27cSBarry Smith   if (ja && a->j != *ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
433d42735a1SBarry Smith   if (oshift) {
43470d8d27cSBarry Smith     if (!ia) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inia[] argument");
43570d8d27cSBarry Smith     if (!ja) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"If oshift then you must passed in inja[] argument");
436d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
437d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
438d42735a1SBarry Smith       (*ja)[i]--;
439d42735a1SBarry Smith     }
440d42735a1SBarry Smith   }
4416cd91f64SBarry Smith   PetscFunctionReturn(0);
4426cd91f64SBarry Smith }
443b97cf49bSBarry Smith 
44417667f90SBarry Smith #undef __FUNCT__
44517667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj"
44619fd82e9SBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
44717667f90SBarry Smith {
44817667f90SBarry Smith   Mat               B;
44917667f90SBarry Smith   PetscErrorCode    ierr;
45017667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
45117667f90SBarry Smith   const PetscInt    *rj;
45217667f90SBarry Smith   const PetscScalar *ra;
45317667f90SBarry Smith   MPI_Comm          comm;
45417667f90SBarry Smith 
45517667f90SBarry Smith   PetscFunctionBegin;
4560298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
4570298fd71SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
4580298fd71SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
45917667f90SBarry Smith 
46017667f90SBarry Smith   /* count the number of nonzeros per row */
46117667f90SBarry Smith   for (i=0; i<m; i++) {
4620298fd71SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
4635ee9ba1cSJed Brown     for (j=0; j<len; j++) {
4645ee9ba1cSJed Brown       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
4655ee9ba1cSJed Brown     }
46617667f90SBarry Smith     nzeros += len;
4670f2063bfSTobin Isaac     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
46817667f90SBarry Smith   }
46917667f90SBarry Smith 
47017667f90SBarry Smith   /* malloc space for nonzeros */
471854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr);
472854ce69bSBarry Smith   ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr);
473854ce69bSBarry Smith   ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr);
47417667f90SBarry Smith 
47517667f90SBarry Smith   nzeros = 0;
47617667f90SBarry Smith   ia[0]  = 0;
47717667f90SBarry Smith   for (i=0; i<m; i++) {
47817667f90SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
47917667f90SBarry Smith     cnt  = 0;
48017667f90SBarry Smith     for (j=0; j<len; j++) {
4815ee9ba1cSJed Brown       if (rj[j] != i+rstart) { /* if not diagonal */
48217667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
48317667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
48417667f90SBarry Smith       }
4855ee9ba1cSJed Brown     }
48617667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
48717667f90SBarry Smith     nzeros += cnt;
48817667f90SBarry Smith     ia[i+1] = nzeros;
48917667f90SBarry Smith   }
49017667f90SBarry Smith 
49117667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
49217667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
49358ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
49417667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
49517667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
49617667f90SBarry Smith 
497511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
49828be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
49917667f90SBarry Smith   } else {
50017667f90SBarry Smith     *newmat = B;
50117667f90SBarry Smith   }
50217667f90SBarry Smith   PetscFunctionReturn(0);
50317667f90SBarry Smith }
50417667f90SBarry Smith 
505b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
5060b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
5073eda8832SBarry Smith                                        MatGetRow_MPIAdj,
5083eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
509b97cf49bSBarry Smith                                        0,
51097304618SKris Buschelman                                 /* 4*/ 0,
511b97cf49bSBarry Smith                                        0,
512b97cf49bSBarry Smith                                        0,
513b97cf49bSBarry Smith                                        0,
5140b3b1487SBarry Smith                                        0,
5150b3b1487SBarry Smith                                        0,
51697304618SKris Buschelman                                 /*10*/ 0,
5170b3b1487SBarry Smith                                        0,
5180b3b1487SBarry Smith                                        0,
5190b3b1487SBarry Smith                                        0,
5200b3b1487SBarry Smith                                        0,
52197304618SKris Buschelman                                 /*15*/ 0,
5223eda8832SBarry Smith                                        MatEqual_MPIAdj,
5230b3b1487SBarry Smith                                        0,
5240b3b1487SBarry Smith                                        0,
5250b3b1487SBarry Smith                                        0,
52697304618SKris Buschelman                                 /*20*/ 0,
5270b3b1487SBarry Smith                                        0,
5283eda8832SBarry Smith                                        MatSetOption_MPIAdj,
5290b3b1487SBarry Smith                                        0,
530d519adbfSMatthew Knepley                                 /*24*/ 0,
5310b3b1487SBarry Smith                                        0,
5320b3b1487SBarry Smith                                        0,
5330b3b1487SBarry Smith                                        0,
5340b3b1487SBarry Smith                                        0,
535d519adbfSMatthew Knepley                                 /*29*/ 0,
5360b3b1487SBarry Smith                                        0,
537273d9f13SBarry Smith                                        0,
538273d9f13SBarry Smith                                        0,
5390b3b1487SBarry Smith                                        0,
540d519adbfSMatthew Knepley                                 /*34*/ 0,
5410b3b1487SBarry Smith                                        0,
5420b3b1487SBarry Smith                                        0,
5430b3b1487SBarry Smith                                        0,
5440b3b1487SBarry Smith                                        0,
545d519adbfSMatthew Knepley                                 /*39*/ 0,
5465495542aSFande Kong                                        MatGetSubMatrices_MPIAdj,
5470b3b1487SBarry Smith                                        0,
5480b3b1487SBarry Smith                                        0,
5490b3b1487SBarry Smith                                        0,
550d519adbfSMatthew Knepley                                 /*44*/ 0,
5510b3b1487SBarry Smith                                        0,
5527d68702bSBarry Smith                                        MatShift_Basic,
5530b3b1487SBarry Smith                                        0,
5540b3b1487SBarry Smith                                        0,
555d519adbfSMatthew Knepley                                 /*49*/ 0,
5566cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
557d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
558b97cf49bSBarry Smith                                        0,
559b97cf49bSBarry Smith                                        0,
560d519adbfSMatthew Knepley                                 /*54*/ 0,
561b97cf49bSBarry Smith                                        0,
5620752156aSBarry Smith                                        0,
5630752156aSBarry Smith                                        0,
5640b3b1487SBarry Smith                                        0,
565d519adbfSMatthew Knepley                                 /*59*/ 0,
566b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
567b9b97703SBarry Smith                                        MatView_MPIAdj,
56817667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
56997304618SKris Buschelman                                        0,
570d519adbfSMatthew Knepley                                 /*64*/ 0,
57197304618SKris Buschelman                                        0,
57297304618SKris Buschelman                                        0,
57397304618SKris Buschelman                                        0,
57497304618SKris Buschelman                                        0,
575d519adbfSMatthew Knepley                                 /*69*/ 0,
57697304618SKris Buschelman                                        0,
57797304618SKris Buschelman                                        0,
57897304618SKris Buschelman                                        0,
57997304618SKris Buschelman                                        0,
580d519adbfSMatthew Knepley                                 /*74*/ 0,
58197304618SKris Buschelman                                        0,
58297304618SKris Buschelman                                        0,
58397304618SKris Buschelman                                        0,
58497304618SKris Buschelman                                        0,
585d519adbfSMatthew Knepley                                 /*79*/ 0,
58697304618SKris Buschelman                                        0,
58797304618SKris Buschelman                                        0,
58897304618SKris Buschelman                                        0,
589865e5f61SKris Buschelman                                        0,
590d519adbfSMatthew Knepley                                 /*84*/ 0,
591865e5f61SKris Buschelman                                        0,
592865e5f61SKris Buschelman                                        0,
593865e5f61SKris Buschelman                                        0,
594865e5f61SKris Buschelman                                        0,
595d519adbfSMatthew Knepley                                 /*89*/ 0,
596865e5f61SKris Buschelman                                        0,
597865e5f61SKris Buschelman                                        0,
598865e5f61SKris Buschelman                                        0,
599865e5f61SKris Buschelman                                        0,
600d519adbfSMatthew Knepley                                 /*94*/ 0,
601865e5f61SKris Buschelman                                        0,
602865e5f61SKris Buschelman                                        0,
6033964eb88SJed Brown                                        0,
6043964eb88SJed Brown                                        0,
6053964eb88SJed Brown                                 /*99*/ 0,
6063964eb88SJed Brown                                        0,
6073964eb88SJed Brown                                        0,
6083964eb88SJed Brown                                        0,
6093964eb88SJed Brown                                        0,
6103964eb88SJed Brown                                /*104*/ 0,
6113964eb88SJed Brown                                        0,
6123964eb88SJed Brown                                        0,
6133964eb88SJed Brown                                        0,
6143964eb88SJed Brown                                        0,
6153964eb88SJed Brown                                /*109*/ 0,
6163964eb88SJed Brown                                        0,
6173964eb88SJed Brown                                        0,
6183964eb88SJed Brown                                        0,
6193964eb88SJed Brown                                        0,
6203964eb88SJed Brown                                /*114*/ 0,
6213964eb88SJed Brown                                        0,
6223964eb88SJed Brown                                        0,
6233964eb88SJed Brown                                        0,
6243964eb88SJed Brown                                        0,
6253964eb88SJed Brown                                /*119*/ 0,
6263964eb88SJed Brown                                        0,
6273964eb88SJed Brown                                        0,
6283964eb88SJed Brown                                        0,
6293964eb88SJed Brown                                        0,
6303964eb88SJed Brown                                /*124*/ 0,
6313964eb88SJed Brown                                        0,
6323964eb88SJed Brown                                        0,
6333964eb88SJed Brown                                        0,
63440475a53SFande Kong                                        MatGetSubMatricesMPI_MPIAdj,
6353964eb88SJed Brown                                /*129*/ 0,
6363964eb88SJed Brown                                        0,
6373964eb88SJed Brown                                        0,
6383964eb88SJed Brown                                        0,
6393964eb88SJed Brown                                        0,
6403964eb88SJed Brown                                /*134*/ 0,
6413964eb88SJed Brown                                        0,
6423964eb88SJed Brown                                        0,
6433964eb88SJed Brown                                        0,
6443964eb88SJed Brown                                        0,
6453964eb88SJed Brown                                /*139*/ 0,
646f9426fe0SMark Adams                                        0,
6473964eb88SJed Brown                                        0
6483964eb88SJed Brown };
649b97cf49bSBarry Smith 
650a23d5eceSKris Buschelman #undef __FUNCT__
651a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
652f7a08781SBarry Smith static PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
653a23d5eceSKris Buschelman {
654a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj*)B->data;
655dfbe8321SBarry Smith   PetscErrorCode ierr;
6562515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
657b24ad042SBarry Smith   PetscInt ii;
658a23d5eceSKris Buschelman #endif
659a23d5eceSKris Buschelman 
660a23d5eceSKris Buschelman   PetscFunctionBegin;
66126283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
66226283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
66358ec18a5SLisandro Dalcin 
6642515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
665e32f2f54SBarry 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]);
666d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
667e02043d6SBarry 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]);
668a23d5eceSKris Buschelman   }
669d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
670e02043d6SBarry 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]);
671a23d5eceSKris Buschelman   }
672a23d5eceSKris Buschelman #endif
67358ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
674a23d5eceSKris Buschelman 
675a23d5eceSKris Buschelman   b->j      = j;
676a23d5eceSKris Buschelman   b->i      = i;
677a23d5eceSKris Buschelman   b->values = values;
678a23d5eceSKris Buschelman 
679d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
680a23d5eceSKris Buschelman   b->diag      = 0;
681a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
682a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
683a23d5eceSKris Buschelman 
684a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
685a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
686a23d5eceSKris Buschelman   PetscFunctionReturn(0);
687a23d5eceSKris Buschelman }
688b97cf49bSBarry Smith 
689841d17a1SFande Kong #undef __FUNCT__
6909a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj"
691f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B)
6929a3665c6SJed Brown {
6939a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
6949a3665c6SJed Brown   PetscErrorCode ierr;
6959a3665c6SJed Brown   const PetscInt *ranges;
6969a3665c6SJed Brown   MPI_Comm       acomm,bcomm;
6979a3665c6SJed Brown   MPI_Group      agroup,bgroup;
6989a3665c6SJed Brown   PetscMPIInt    i,rank,size,nranks,*ranks;
6999a3665c6SJed Brown 
7009a3665c6SJed Brown   PetscFunctionBegin;
7010298fd71SBarry Smith   *B    = NULL;
702ce94432eSBarry Smith   ierr  = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr);
7039a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&size);CHKERRQ(ierr);
7049a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr);
7059a3665c6SJed Brown   ierr  = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
7069a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
7079a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) nranks++;
7089a3665c6SJed Brown   }
7099a3665c6SJed Brown   if (nranks == size) {         /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
7109a3665c6SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
7119a3665c6SJed Brown     *B   = A;
7129a3665c6SJed Brown     PetscFunctionReturn(0);
7139a3665c6SJed Brown   }
7149a3665c6SJed Brown 
715785e854fSJed Brown   ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr);
7169a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
7179a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i;
7189a3665c6SJed Brown   }
7199a3665c6SJed Brown   ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr);
7209a3665c6SJed Brown   ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr);
7219a3665c6SJed Brown   ierr = PetscFree(ranks);CHKERRQ(ierr);
7229a3665c6SJed Brown   ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr);
7239a3665c6SJed Brown   ierr = MPI_Group_free(&agroup);CHKERRQ(ierr);
7249a3665c6SJed Brown   ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr);
7259a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
7269a3665c6SJed Brown     PetscInt   m,N;
7279a3665c6SJed Brown     Mat_MPIAdj *b;
7280298fd71SBarry Smith     ierr       = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
7290298fd71SBarry Smith     ierr       = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
7309a3665c6SJed Brown     ierr       = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr);
7319a3665c6SJed Brown     b          = (Mat_MPIAdj*)(*B)->data;
7329a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
7339a3665c6SJed Brown     ierr       = MPI_Comm_free(&bcomm);CHKERRQ(ierr);
7349a3665c6SJed Brown   }
7359a3665c6SJed Brown   PetscFunctionReturn(0);
7369a3665c6SJed Brown }
7379a3665c6SJed Brown 
7389a3665c6SJed Brown #undef __FUNCT__
7399a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat"
7409a3665c6SJed Brown /*@
7419a3665c6SJed Brown    MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows
7429a3665c6SJed Brown 
7439a3665c6SJed Brown    Collective
7449a3665c6SJed Brown 
7459a3665c6SJed Brown    Input Arguments:
7469a3665c6SJed Brown .  A - original MPIAdj matrix
7479a3665c6SJed Brown 
7489a3665c6SJed Brown    Output Arguments:
7490298fd71SBarry Smith .  B - matrix on subcommunicator, NULL on ranks that owned zero rows of A
7509a3665c6SJed Brown 
7519a3665c6SJed Brown    Level: developer
7529a3665c6SJed Brown 
7539a3665c6SJed Brown    Note:
7549a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
7559a3665c6SJed Brown 
7569a3665c6SJed Brown    The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed.
7579a3665c6SJed Brown 
7589a3665c6SJed Brown .seealso: MatCreateMPIAdj()
7599a3665c6SJed Brown @*/
7609a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B)
7619a3665c6SJed Brown {
7629a3665c6SJed Brown   PetscErrorCode ierr;
7639a3665c6SJed Brown 
7649a3665c6SJed Brown   PetscFunctionBegin;
7659a3665c6SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
7669a3665c6SJed Brown   ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
7679a3665c6SJed Brown   PetscFunctionReturn(0);
7689a3665c6SJed Brown }
7699a3665c6SJed Brown 
7700bad9183SKris Buschelman /*MC
771fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
7720bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
7730bad9183SKris Buschelman 
7740bad9183SKris Buschelman   Level: beginner
7750bad9183SKris Buschelman 
7760bad9183SKris Buschelman .seealso: MatCreateMPIAdj
7770bad9183SKris Buschelman M*/
7780bad9183SKris Buschelman 
7794a2ae208SSatish Balay #undef __FUNCT__
7804a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
7818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
782273d9f13SBarry Smith {
783273d9f13SBarry Smith   Mat_MPIAdj     *b;
7846849ba73SBarry Smith   PetscErrorCode ierr;
785273d9f13SBarry Smith 
786273d9f13SBarry Smith   PetscFunctionBegin;
787b00a9115SJed Brown   ierr         = PetscNewLog(B,&b);CHKERRQ(ierr);
788b0a32e0cSBarry Smith   B->data      = (void*)b;
789273d9f13SBarry Smith   ierr         = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
790273d9f13SBarry Smith   B->assembled = PETSC_FALSE;
791273d9f13SBarry Smith 
792bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
793bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr);
79417667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
795273d9f13SBarry Smith   PetscFunctionReturn(0);
796273d9f13SBarry Smith }
797273d9f13SBarry Smith 
7984a2ae208SSatish Balay #undef __FUNCT__
7994a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
800a23d5eceSKris Buschelman /*@C
801a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
802a23d5eceSKris Buschelman 
8033f9fe445SBarry Smith    Logically Collective on MPI_Comm
804a23d5eceSKris Buschelman 
805a23d5eceSKris Buschelman    Input Parameters:
806a23d5eceSKris Buschelman +  A - the matrix
807a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
808a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
809a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
810a23d5eceSKris Buschelman -  values - [optional] edge weights
811a23d5eceSKris Buschelman 
812a23d5eceSKris Buschelman    Level: intermediate
813a23d5eceSKris Buschelman 
814a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
815a23d5eceSKris Buschelman @*/
8167087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
817273d9f13SBarry Smith {
8184ac538c5SBarry Smith   PetscErrorCode ierr;
819273d9f13SBarry Smith 
820273d9f13SBarry Smith   PetscFunctionBegin;
8214ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
822273d9f13SBarry Smith   PetscFunctionReturn(0);
823273d9f13SBarry Smith }
824273d9f13SBarry Smith 
8254a2ae208SSatish Balay #undef __FUNCT__
8264a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
827b97cf49bSBarry Smith /*@C
8283eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
829b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
830b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
831b97cf49bSBarry Smith 
832ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
833ef5ee4d1SLois Curfman McInnes 
834b97cf49bSBarry Smith    Input Parameters:
835c2e958feSBarry Smith +  comm - MPI communicator
8360752156aSBarry Smith .  m - number of local rows
83758ec18a5SLisandro Dalcin .  N - number of global columns
838b97cf49bSBarry Smith .  i - the indices into j for the start of each row
839329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
840ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
841329f5518SBarry Smith -  values -[optional] edge weights
842b97cf49bSBarry Smith 
843b97cf49bSBarry Smith    Output Parameter:
844b97cf49bSBarry Smith .  A - the matrix
845b97cf49bSBarry Smith 
84615091d37SBarry Smith    Level: intermediate
84715091d37SBarry Smith 
8484bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
8494bc6d8c0SBarry Smith    MatSetValues().
850c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
851ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
8521198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
853327e1a73SBarry Smith    Should not include the matrix diagonals.
854b97cf49bSBarry Smith 
855e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
856ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
857ededeb1aSvictorle 
858ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
859b97cf49bSBarry Smith 
860e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
861b97cf49bSBarry Smith @*/
8627087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
863b97cf49bSBarry Smith {
864dfbe8321SBarry Smith   PetscErrorCode ierr;
865b97cf49bSBarry Smith 
866433994e6SBarry Smith   PetscFunctionBegin;
867f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
86858ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
869273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
870273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
8713a40ed3dSBarry Smith   PetscFunctionReturn(0);
872b97cf49bSBarry Smith }
873b97cf49bSBarry Smith 
874c2e958feSBarry Smith 
875433994e6SBarry Smith 
876433994e6SBarry Smith 
877433994e6SBarry Smith 
878433994e6SBarry Smith 
879433994e6SBarry Smith 
880