xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 30de9b251c71d073cc9d88242eefe8db4dcc3e99)
1 /*$Id: mpiadj.c,v 1.66 2001/08/07 03:02:59 balay Exp $*/
2 
3 /*
4     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
5 */
6 #include "src/mat/impls/adj/mpi/mpiadj.h"
7 #include "petscsys.h"
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "MatView_MPIAdj_ASCII"
11 int MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
12 {
13   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
14   int               ierr,i,j,m = A->m;
15   char              *name;
16   PetscViewerFormat format;
17 
18   PetscFunctionBegin;
19   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
20   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
21   if (format == PETSC_VIEWER_ASCII_INFO) {
22     PetscFunctionReturn(0);
23   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
24     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
25   } else {
26     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
27     for (i=0; i<m; i++) {
28       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr);
29       for (j=a->i[i]; j<a->i[i+1]; j++) {
30         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr);
31       }
32       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
33     }
34     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
35   }
36   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "MatView_MPIAdj"
42 int MatView_MPIAdj(Mat A,PetscViewer viewer)
43 {
44   int        ierr;
45   PetscTruth isascii;
46 
47   PetscFunctionBegin;
48   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
49   if (isascii) {
50     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
51   } else {
52     SETERRQ1(1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name);
53   }
54   PetscFunctionReturn(0);
55 }
56 
57 #undef __FUNCT__
58 #define __FUNCT__ "MatDestroy_MPIAdj"
59 int MatDestroy_MPIAdj(Mat mat)
60 {
61   Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data;
62   int        ierr;
63 
64   PetscFunctionBegin;
65 #if defined(PETSC_USE_LOG)
66   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz);
67 #endif
68   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
69   if (a->freeaij) {
70     ierr = PetscFree(a->i);CHKERRQ(ierr);
71     ierr = PetscFree(a->j);CHKERRQ(ierr);
72     if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);}
73   }
74   ierr = PetscFree(a->rowners);CHKERRQ(ierr);
75   ierr = PetscFree(a);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "MatSetOption_MPIAdj"
81 int MatSetOption_MPIAdj(Mat A,MatOption op)
82 {
83   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
84 
85   PetscFunctionBegin;
86   switch (op) {
87   case MAT_STRUCTURALLY_SYMMETRIC:
88     a->symmetric = PETSC_TRUE;
89     break;
90   default:
91     PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
92     break;
93   }
94   PetscFunctionReturn(0);
95 }
96 
97 
98 /*
99      Adds diagonal pointers to sparse matrix structure.
100 */
101 
102 #undef __FUNCT__
103 #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
104 int MatMarkDiagonal_MPIAdj(Mat A)
105 {
106   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
107   int        i,j,*diag,m = A->m,ierr;
108 
109   PetscFunctionBegin;
110   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
111   PetscLogObjectMemory(A,(m+1)*sizeof(int));
112   for (i=0; i<A->m; i++) {
113     for (j=a->i[i]; j<a->i[i+1]; j++) {
114       if (a->j[j] == i) {
115         diag[i] = j;
116         break;
117       }
118     }
119   }
120   a->diag = diag;
121   PetscFunctionReturn(0);
122 }
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "MatGetRow_MPIAdj"
126 int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v)
127 {
128   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
129   int        *itmp;
130 
131   PetscFunctionBegin;
132   row -= a->rstart;
133 
134   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
135 
136   *nz = a->i[row+1] - a->i[row];
137   if (v) *v = PETSC_NULL;
138   if (idx) {
139     itmp = a->j + a->i[row];
140     if (*nz) {
141       *idx = itmp;
142     }
143     else *idx = 0;
144   }
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "MatRestoreRow_MPIAdj"
150 int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v)
151 {
152   PetscFunctionBegin;
153   PetscFunctionReturn(0);
154 }
155 
156 #undef __FUNCT__
157 #define __FUNCT__ "MatGetBlockSize_MPIAdj"
158 int MatGetBlockSize_MPIAdj(Mat A,int *bs)
159 {
160   PetscFunctionBegin;
161   *bs = 1;
162   PetscFunctionReturn(0);
163 }
164 
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "MatEqual_MPIAdj"
168 int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
169 {
170   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
171   int         ierr;
172   PetscTruth  flag;
173 
174   PetscFunctionBegin;
175   ierr = PetscTypeCompare((PetscObject)B,MATMPIADJ,&flag);CHKERRQ(ierr);
176   if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");
177 
178   /* If the  matrix dimensions are not equal,or no of nonzeros */
179   if ((A->m != B->m) ||(a->nz != b->nz)) {
180     flag = PETSC_FALSE;
181   }
182 
183   /* if the a->i are the same */
184   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr);
185 
186   /* if a->j are the same */
187   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr);
188 
189   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
190   PetscFunctionReturn(0);
191 }
192 
193 #undef __FUNCT__
194 #define __FUNCT__ "MatGetRowIJ_MPIAdj"
195 int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
196 {
197   int        ierr,size,i;
198   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
199 
200   PetscFunctionBegin;
201   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
202   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
203   *m    = A->m;
204   *ia   = a->i;
205   *ja   = a->j;
206   *done = PETSC_TRUE;
207   if (oshift) {
208     for (i=0; i<(*ia)[*m]; i++) {
209       (*ja)[i]++;
210     }
211     for (i=0; i<=(*m); i++) (*ia)[i]++;
212   }
213   PetscFunctionReturn(0);
214 }
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
218 int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
219 {
220   int        i;
221   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
222 
223   PetscFunctionBegin;
224   if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()");
225   if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()");
226   if (oshift) {
227     for (i=0; i<=(*m); i++) (*ia)[i]--;
228     for (i=0; i<(*ia)[*m]; i++) {
229       (*ja)[i]--;
230     }
231   }
232   PetscFunctionReturn(0);
233 }
234 
235 /* -------------------------------------------------------------------*/
236 static struct _MatOps MatOps_Values = {0,
237        MatGetRow_MPIAdj,
238        MatRestoreRow_MPIAdj,
239        0,
240        0,
241        0,
242        0,
243        0,
244        0,
245        0,
246        0,
247        0,
248        0,
249        0,
250        0,
251        0,
252        MatEqual_MPIAdj,
253        0,
254        0,
255        0,
256        0,
257        0,
258        0,
259        MatSetOption_MPIAdj,
260        0,
261        0,
262        0,
263        0,
264        0,
265        0,
266        0,
267        0,
268        0,
269        0,
270        0,
271        0,
272        0,
273        0,
274        0,
275        0,
276        0,
277        0,
278        0,
279        0,
280        0,
281        0,
282        0,
283        0,
284        0,
285        0,
286        MatGetBlockSize_MPIAdj,
287        MatGetRowIJ_MPIAdj,
288        MatRestoreRowIJ_MPIAdj,
289        0,
290        0,
291        0,
292        0,
293        0,
294        0,
295        0,
296        0,
297        MatDestroy_MPIAdj,
298        MatView_MPIAdj,
299        MatGetPetscMaps_Petsc};
300 
301 
302 EXTERN_C_BEGIN
303 #undef __FUNCT__
304 #define __FUNCT__ "MatCreate_MPIAdj"
305 int MatCreate_MPIAdj(Mat B)
306 {
307   Mat_MPIAdj *b;
308   int        ii,ierr,size,rank;
309 
310   PetscFunctionBegin;
311   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
312   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
313 
314   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
315   B->data             = (void*)b;
316   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
317   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
318   B->factor           = 0;
319   B->lupivotthreshold = 1.0;
320   B->mapping          = 0;
321   B->assembled        = PETSC_FALSE;
322 
323   ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr);
324   B->n = B->N;
325 
326   /* the information in the maps duplicates the information computed below, eventually
327      we should remove the duplicate information that is not contained in the maps */
328   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
329   /* we don't know the "local columns" so just use the row information :-(*/
330   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
331 
332   ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr);
333   PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
334   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
335   b->rowners[0] = 0;
336   for (ii=2; ii<=size; ii++) {
337     b->rowners[ii] += b->rowners[ii-1];
338   }
339   b->rstart = b->rowners[rank];
340   b->rend   = b->rowners[rank+1];
341 
342   PetscFunctionReturn(0);
343 }
344 EXTERN_C_END
345 
346 #undef __FUNCT__
347 #define __FUNCT__ "MatMPIAdjSetPreallocation"
348 int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values)
349 {
350   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
351   int        ierr;
352 #if defined(PETSC_USE_BOPT_g)
353   int        ii;
354 #endif
355 
356   PetscFunctionBegin;
357   B->preallocated = PETSC_TRUE;
358 #if defined(PETSC_USE_BOPT_g)
359   if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]);
360   for (ii=1; ii<B->m; ii++) {
361     if (i[ii] < 0 || i[ii] < i[ii-1]) {
362       SETERRQ4(1,"i[%d]=%d index is out of range: i[%d]=%d",ii,i[ii],ii-1,i[ii-1]);
363     }
364   }
365   for (ii=0; ii<i[B->m]; ii++) {
366     if (j[ii] < 0 || j[ii] >= B->N) {
367       SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]);
368     }
369   }
370 #endif
371 
372   b->j      = j;
373   b->i      = i;
374   b->values = values;
375 
376   b->nz               = i[B->m];
377   b->diag             = 0;
378   b->symmetric        = PETSC_FALSE;
379   b->freeaij          = PETSC_TRUE;
380 
381   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
382   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
383   PetscFunctionReturn(0);
384 }
385 
386 #undef __FUNCT__
387 #define __FUNCT__ "MatCreateMPIAdj"
388 /*@C
389    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
390    The matrix does not have numerical values associated with it, but is
391    intended for ordering (to reduce bandwidth etc) and partitioning.
392 
393    Collective on MPI_Comm
394 
395    Input Parameters:
396 +  comm - MPI communicator
397 .  m - number of local rows
398 .  n - number of columns
399 .  i - the indices into j for the start of each row
400 .  j - the column indices for each row (sorted for each row).
401        The indices in i and j start with zero (NOT with one).
402 -  values -[optional] edge weights
403 
404    Output Parameter:
405 .  A - the matrix
406 
407    Level: intermediate
408 
409    Notes: This matrix object does not support most matrix operations, include
410    MatSetValues().
411    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
412    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
413     call from Fortran you need not create the arrays with PetscMalloc().
414    Should not include the matrix diagonals.
415 
416    If you already have a matrix, you can create its adjacency matrix by a call
417    to MatConvert, specifying a type of MATMPIADJ.
418 
419    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
420 
421 .seealso: MatCreate(), MatConvert(), MatGetOrdering()
422 @*/
423 int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A)
424 {
425   int        ierr;
426 
427   PetscFunctionBegin;
428   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
429   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
430   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
431   PetscFunctionReturn(0);
432 }
433 
434 EXTERN_C_BEGIN
435 #undef __FUNCT__
436 #define __FUNCT__ "MatConvertTo_MPIAdj"
437 int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B)
438 {
439   int          i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a;
440   PetscScalar  *ra;
441   MPI_Comm     comm;
442 
443   PetscFunctionBegin;
444   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
445   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
446   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
447 
448   /* count the number of nonzeros per row */
449   for (i=0; i<m; i++) {
450     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
451     for (j=0; j<len; j++) {
452       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
453     }
454     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
455     nzeros += len;
456   }
457 
458   /* malloc space for nonzeros */
459   ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr);
460   ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr);
461   ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr);
462 
463   nzeros = 0;
464   ia[0]  = 0;
465   for (i=0; i<m; i++) {
466     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
467     cnt     = 0;
468     for (j=0; j<len; j++) {
469       if (rj[j] != i+rstart) { /* if not diagonal */
470         a[nzeros+cnt]    = (int) PetscAbsScalar(ra[j]);
471         ja[nzeros+cnt++] = rj[j];
472       }
473     }
474     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
475     nzeros += cnt;
476     ia[i+1] = nzeros;
477   }
478 
479   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
480   ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr);
481 
482   PetscFunctionReturn(0);
483 }
484 EXTERN_C_END
485 
486 
487 
488 
489 
490