xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 5cd905551efc76d042fcd2cc746dd1c50567705b)
1b97cf49bSBarry Smith #ifdef PETSC_RCS_HEADER
2*5cd90555SBarry Smith static char vcid[] = "$Id: mpiadj.c,v 1.9 1998/04/03 23:15:52 bsmith Exp bsmith $";
3b97cf49bSBarry Smith #endif
4b97cf49bSBarry Smith 
5b97cf49bSBarry Smith /*
6b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
7b97cf49bSBarry Smith */
8b97cf49bSBarry Smith 
9b97cf49bSBarry Smith #include "pinclude/pviewer.h"
10b97cf49bSBarry Smith #include "sys.h"
110752156aSBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h"
12b97cf49bSBarry Smith 
13b97cf49bSBarry Smith 
14b97cf49bSBarry Smith #undef __FUNC__
150752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj_ASCII"
160752156aSBarry Smith extern int MatView_MPIAdj_ASCII(Mat A,Viewer viewer)
17b97cf49bSBarry Smith {
180752156aSBarry Smith   Mat_MPIAdj  *a = (Mat_MPIAdj *) A->data;
19b97cf49bSBarry Smith   int         ierr, i,j, m = a->m,  format;
20b97cf49bSBarry Smith   FILE        *fd;
21b97cf49bSBarry Smith   char        *outputname;
220752156aSBarry Smith   MPI_Comm    comm = A->comm;
23b97cf49bSBarry Smith 
24b97cf49bSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
25b97cf49bSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
26b97cf49bSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
27b97cf49bSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO) {
283a40ed3dSBarry Smith     PetscFunctionReturn(0);
290752156aSBarry Smith   } else {
30b97cf49bSBarry Smith     for ( i=0; i<m; i++ ) {
310752156aSBarry Smith       PetscSynchronizedFPrintf(comm,fd,"row %d:",i+a->rstart);
32b97cf49bSBarry Smith       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
330752156aSBarry Smith         PetscSynchronizedFPrintf(comm,fd," %d ",a->j[j]);
340752156aSBarry Smith       }
350752156aSBarry Smith       PetscSynchronizedFPrintf(comm,fd,"\n");
36b97cf49bSBarry Smith     }
37b97cf49bSBarry Smith   }
380752156aSBarry Smith   PetscSynchronizedFlush(comm);
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
40b97cf49bSBarry Smith }
41b97cf49bSBarry Smith 
42b97cf49bSBarry Smith #undef __FUNC__
430752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj"
44e1311b90SBarry Smith int MatView_MPIAdj(Mat A,Viewer viewer)
45b97cf49bSBarry Smith {
46b97cf49bSBarry Smith   ViewerType  vtype;
47b97cf49bSBarry Smith   int         ierr;
48b97cf49bSBarry Smith 
49b97cf49bSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
500752156aSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
513a40ed3dSBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
52*5cd90555SBarry Smith   } else {
53*5cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
54b97cf49bSBarry Smith   }
553a40ed3dSBarry Smith   PetscFunctionReturn(0);
56b97cf49bSBarry Smith }
57b97cf49bSBarry Smith 
58b97cf49bSBarry Smith #undef __FUNC__
590752156aSBarry Smith #define __FUNC__ "MatDestroy_MPIAdj"
60e1311b90SBarry Smith int MatDestroy_MPIAdj(Mat A)
61b97cf49bSBarry Smith {
620752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
63b97cf49bSBarry Smith 
643a40ed3dSBarry Smith #if defined(USE_PETSC_LOG)
65e1311b90SBarry Smith   PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz);
66b97cf49bSBarry Smith #endif
67b97cf49bSBarry Smith   if (a->diag) PetscFree(a->diag);
68b97cf49bSBarry Smith   PetscFree(a->i);
69b97cf49bSBarry Smith   PetscFree(a->j);
700752156aSBarry Smith   PetscFree(a->rowners);
71b97cf49bSBarry Smith   PetscFree(a);
72b97cf49bSBarry Smith 
73b97cf49bSBarry Smith   PLogObjectDestroy(A);
74b97cf49bSBarry Smith   PetscHeaderDestroy(A);
753a40ed3dSBarry Smith   PetscFunctionReturn(0);
76b97cf49bSBarry Smith }
77b97cf49bSBarry Smith 
78b97cf49bSBarry Smith 
79b97cf49bSBarry Smith #undef __FUNC__
800752156aSBarry Smith #define __FUNC__ "MatSetOption_MPIAdj"
810752156aSBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op)
82b97cf49bSBarry Smith {
830752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
84b97cf49bSBarry Smith 
85b97cf49bSBarry Smith   if (op == MAT_STRUCTURALLY_SYMMETRIC) {
86b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
87b97cf49bSBarry Smith   } else {
88981c4779SBarry Smith     PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
89b97cf49bSBarry Smith   }
903a40ed3dSBarry Smith   PetscFunctionReturn(0);
91b97cf49bSBarry Smith }
92b97cf49bSBarry Smith 
93b97cf49bSBarry Smith 
94b97cf49bSBarry Smith /*
95b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
96b97cf49bSBarry Smith */
97b97cf49bSBarry Smith 
98b97cf49bSBarry Smith #undef __FUNC__
990752156aSBarry Smith #define __FUNC__ "MatMarkDiag_MPIAdj"
1000752156aSBarry Smith int MatMarkDiag_MPIAdj(Mat A)
101b97cf49bSBarry Smith {
1020752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
103b97cf49bSBarry Smith   int        i,j, *diag, m = a->m;
104b97cf49bSBarry Smith 
105b97cf49bSBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
106b97cf49bSBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
107b97cf49bSBarry Smith   for ( i=0; i<a->m; i++ ) {
108b97cf49bSBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
109b97cf49bSBarry Smith       if (a->j[j] == i) {
110b97cf49bSBarry Smith         diag[i] = j;
111b97cf49bSBarry Smith         break;
112b97cf49bSBarry Smith       }
113b97cf49bSBarry Smith     }
114b97cf49bSBarry Smith   }
115b97cf49bSBarry Smith   a->diag = diag;
1163a40ed3dSBarry Smith   PetscFunctionReturn(0);
117b97cf49bSBarry Smith }
118b97cf49bSBarry Smith 
119b97cf49bSBarry Smith #undef __FUNC__
1200752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj"
1210752156aSBarry Smith int MatGetSize_MPIAdj(Mat A,int *m,int *n)
122b97cf49bSBarry Smith {
1230752156aSBarry Smith   if (m) *m = A->M;
1240752156aSBarry Smith   if (n) *n = A->N;
1253a40ed3dSBarry Smith   PetscFunctionReturn(0);
126b97cf49bSBarry Smith }
127b97cf49bSBarry Smith 
128b97cf49bSBarry Smith #undef __FUNC__
1290752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj"
1300752156aSBarry Smith int MatGetLocalSize_MPIAdj(Mat A,int *m,int *n)
131b97cf49bSBarry Smith {
1320752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
1330752156aSBarry Smith   if (m) *m = a->m;
1340752156aSBarry Smith   if (n) *n = A->N;
1353a40ed3dSBarry Smith   PetscFunctionReturn(0);
136b97cf49bSBarry Smith }
137b97cf49bSBarry Smith 
138b97cf49bSBarry Smith #undef __FUNC__
1390752156aSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAdj"
1400752156aSBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n)
141b97cf49bSBarry Smith {
1420752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
1430752156aSBarry Smith   *m = a->rstart; *n = a->rend;
1443a40ed3dSBarry Smith   PetscFunctionReturn(0);
145b97cf49bSBarry Smith }
146b97cf49bSBarry Smith 
147b97cf49bSBarry Smith #undef __FUNC__
1480752156aSBarry Smith #define __FUNC__ "MatGetRow_MPIAdj"
1490752156aSBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
150b97cf49bSBarry Smith {
1510752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *) A->data;
152b97cf49bSBarry Smith   int        *itmp;
153b97cf49bSBarry Smith 
1540752156aSBarry Smith   row -= a->rstart;
1550752156aSBarry Smith 
156a8c6a408SBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
157b97cf49bSBarry Smith 
158b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
159b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
160b97cf49bSBarry Smith   if (idx) {
161b97cf49bSBarry Smith     itmp = a->j + a->i[row];
162b97cf49bSBarry Smith     if (*nz) {
163b97cf49bSBarry Smith       *idx = itmp;
164b97cf49bSBarry Smith     }
165b97cf49bSBarry Smith     else *idx = 0;
166b97cf49bSBarry Smith   }
1673a40ed3dSBarry Smith   PetscFunctionReturn(0);
168b97cf49bSBarry Smith }
169b97cf49bSBarry Smith 
170b97cf49bSBarry Smith #undef __FUNC__
1710752156aSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAdj"
1720752156aSBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v)
173b97cf49bSBarry Smith {
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
175b97cf49bSBarry Smith }
176b97cf49bSBarry Smith 
177b97cf49bSBarry Smith #undef __FUNC__
1780752156aSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAdj"
1790752156aSBarry Smith int MatGetBlockSize_MPIAdj(Mat A, int *bs)
180b97cf49bSBarry Smith {
181b97cf49bSBarry Smith   *bs = 1;
1823a40ed3dSBarry Smith   PetscFunctionReturn(0);
183b97cf49bSBarry Smith }
184b97cf49bSBarry Smith 
185b97cf49bSBarry Smith 
186b97cf49bSBarry Smith #undef __FUNC__
1870752156aSBarry Smith #define __FUNC__ "MatEqual_MPIAdj"
1880752156aSBarry Smith int MatEqual_MPIAdj(Mat A,Mat B, PetscTruth* flg)
189b97cf49bSBarry Smith {
1900752156aSBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data;
191ca161407SBarry Smith  int         flag = 1,ierr;
192b97cf49bSBarry Smith 
193a8c6a408SBarry Smith   if (B->type != MATMPIADJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
194b97cf49bSBarry Smith 
195b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal, or no of nonzeros */
1960752156aSBarry Smith   if ((a->m != b->m ) ||( a->nz != b->nz)) {
1970752156aSBarry Smith     flag = 0;
198b97cf49bSBarry Smith   }
199b97cf49bSBarry Smith 
200b97cf49bSBarry Smith   /* if the a->i are the same */
201b97cf49bSBarry Smith   if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) {
2020752156aSBarry Smith     flag = 0;
203b97cf49bSBarry Smith   }
204b97cf49bSBarry Smith 
205b97cf49bSBarry Smith   /* if a->j are the same */
206b97cf49bSBarry Smith   if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) {
2070752156aSBarry Smith     flag = 0;
208b97cf49bSBarry Smith   }
209b97cf49bSBarry Smith 
210ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
2110752156aSBarry Smith 
2120752156aSBarry Smith 
2133a40ed3dSBarry Smith   PetscFunctionReturn(0);
214b97cf49bSBarry Smith }
215b97cf49bSBarry Smith 
216b97cf49bSBarry Smith 
217b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
218b97cf49bSBarry Smith static struct _MatOps MatOps = {0,
2190752156aSBarry Smith        MatGetRow_MPIAdj,MatRestoreRow_MPIAdj,
220b97cf49bSBarry Smith        0,0,
221b97cf49bSBarry Smith        0,0,
222b97cf49bSBarry Smith        0,0,
223b97cf49bSBarry Smith        0,0,
224b97cf49bSBarry Smith        0,0,
225b97cf49bSBarry Smith        0,
226b97cf49bSBarry Smith        0,
2270752156aSBarry Smith        0,MatEqual_MPIAdj,
228b97cf49bSBarry Smith        0,0,0,
229b97cf49bSBarry Smith        0,0,
230b97cf49bSBarry Smith        0,
2310752156aSBarry Smith        MatSetOption_MPIAdj,0,0,
232b97cf49bSBarry Smith        0,0,0,0,
2330752156aSBarry Smith        MatGetSize_MPIAdj,MatGetLocalSize_MPIAdj,MatGetOwnershipRange_MPIAdj,
234b97cf49bSBarry Smith        0,0,
235b97cf49bSBarry Smith        0,0,
236b97cf49bSBarry Smith        0,0,0,
237b97cf49bSBarry Smith        0,0,0,
2380752156aSBarry Smith        0,0,
239b97cf49bSBarry Smith        0,0,
240b97cf49bSBarry Smith        0,
241b97cf49bSBarry Smith        0,0,0,
242b97cf49bSBarry Smith        0,
2430752156aSBarry Smith        MatGetBlockSize_MPIAdj,
244b97cf49bSBarry Smith        0,
245b97cf49bSBarry Smith        0,
246b97cf49bSBarry Smith        0,
247b97cf49bSBarry Smith        0,
248b97cf49bSBarry Smith        0,
2490752156aSBarry Smith        0,
2500752156aSBarry Smith        0,
2510752156aSBarry Smith        0};
252b97cf49bSBarry Smith 
253b97cf49bSBarry Smith 
254b97cf49bSBarry Smith #undef __FUNC__
2550752156aSBarry Smith #define __FUNC__ "MatCreateMPIAdj"
256b97cf49bSBarry Smith /*@C
2570752156aSBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
258b97cf49bSBarry Smith      The matrix does not have numerical values associated with it, but is
259b97cf49bSBarry Smith      intended for ordering (to reduce bandwidth etc) and partitioning.
260b97cf49bSBarry Smith 
261b97cf49bSBarry Smith    Input Parameters:
262b97cf49bSBarry Smith .  comm - MPI communicator, set to PETSC_COMM_SELF
2630752156aSBarry Smith .  m - number of local rows
264b97cf49bSBarry Smith .  n - number of columns
265b97cf49bSBarry Smith .  i - the indices into j for the start of each row
266b97cf49bSBarry Smith .  j - the column indices for each row (sorted for each row)
267b97cf49bSBarry Smith        The indices in i and j start with zero NOT one.
268b97cf49bSBarry Smith 
269b97cf49bSBarry Smith    Output Parameter:
270b97cf49bSBarry Smith .  A - the matrix
271b97cf49bSBarry Smith 
2720752156aSBarry Smith    Notes: You must NOT free the ii and jj arrays yourself. PETSc will free them
273b97cf49bSBarry Smith    when the matrix is destroyed.
274b97cf49bSBarry Smith 
275b58e66f3SSatish Balay .  MatSetOption() possible values - MAT_STRUCTURALLY_SYMMETRIC
276b97cf49bSBarry Smith 
2770752156aSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetReordering()
278b97cf49bSBarry Smith @*/
2790752156aSBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j, Mat *A)
280b97cf49bSBarry Smith {
281b97cf49bSBarry Smith   Mat        B;
2820752156aSBarry Smith   Mat_MPIAdj *b;
2830752156aSBarry Smith   int        ii,ierr, flg,size,rank;
284b97cf49bSBarry Smith 
285b97cf49bSBarry Smith   MPI_Comm_size(comm,&size);
2860752156aSBarry Smith   MPI_Comm_rank(comm,&rank);
287b97cf49bSBarry Smith 
288b97cf49bSBarry Smith   *A                  = 0;
289f830108cSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIADJ,comm,MatDestroy,MatView);
290b97cf49bSBarry Smith   PLogObjectCreate(B);
2910752156aSBarry Smith   B->data             = (void *) (b = PetscNew(Mat_MPIAdj)); CHKPTRQ(b);
2920752156aSBarry Smith   PetscMemzero(b,sizeof(Mat_MPIAdj));
293f830108cSBarry Smith   PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps));
294e1311b90SBarry Smith   B->ops->destroy          = MatDestroy_MPIAdj;
295e1311b90SBarry Smith   B->ops->view             = MatView_MPIAdj;
296b97cf49bSBarry Smith   B->factor           = 0;
297b97cf49bSBarry Smith   B->lupivotthreshold = 1.0;
298b97cf49bSBarry Smith   B->mapping          = 0;
299b97cf49bSBarry Smith   B->assembled        = PETSC_FALSE;
300b97cf49bSBarry Smith 
3010752156aSBarry Smith   b->m = m; B->m = m;
302ca161407SBarry Smith   ierr = MPI_Allreduce(&m,&B->M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
3030752156aSBarry Smith   B->n = n; B->N = n;
3040752156aSBarry Smith 
3050752156aSBarry Smith   b->rowners = (int *) PetscMalloc((size+1)*sizeof(int)); CHKPTRQ(b->rowners);
3060752156aSBarry Smith   PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
307ca161407SBarry Smith   ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
3080752156aSBarry Smith   b->rowners[0] = 0;
3090752156aSBarry Smith   for ( ii=2; ii<=size; ii++ ) {
3100752156aSBarry Smith     b->rowners[ii] += b->rowners[ii-1];
3110752156aSBarry Smith   }
3120752156aSBarry Smith   b->rstart = b->rowners[rank];
3130752156aSBarry Smith   b->rend   = b->rowners[rank+1];
314b97cf49bSBarry Smith 
315b97cf49bSBarry Smith   b->j  = j;
316b97cf49bSBarry Smith   b->i  = i;
317b97cf49bSBarry Smith 
318b97cf49bSBarry Smith   b->nz               = i[m];
319b97cf49bSBarry Smith   b->diag             = 0;
320b97cf49bSBarry Smith   b->symmetric        = PETSC_FALSE;
321b97cf49bSBarry Smith 
322b97cf49bSBarry Smith   *A = B;
323b97cf49bSBarry Smith 
324b97cf49bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
325b97cf49bSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
326b97cf49bSBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
327b97cf49bSBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3283a40ed3dSBarry Smith   PetscFunctionReturn(0);
329b97cf49bSBarry Smith }
330b97cf49bSBarry Smith 
331b97cf49bSBarry Smith 
332b97cf49bSBarry Smith 
333