xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1b97cf49bSBarry Smith /*
2b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
3b97cf49bSBarry Smith */
43eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h"
5325e03aeSBarry Smith #include "petscsys.h"
6b97cf49bSBarry Smith 
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII"
9dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
10b97cf49bSBarry Smith {
113eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
12dfbe8321SBarry Smith   PetscErrorCode ierr;
13dfbe8321SBarry Smith   int i,j,m = A->m;
14fb9695e5SSatish Balay   char              *name;
15ce1411ecSBarry Smith   PetscViewerFormat format;
16b97cf49bSBarry Smith 
17433994e6SBarry Smith   PetscFunctionBegin;
183a7fca6bSBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
19b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
20fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
213a40ed3dSBarry Smith     PetscFunctionReturn(0);
22fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
23ffac6cdbSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
240752156aSBarry Smith   } else {
25b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
26b97cf49bSBarry Smith     for (i=0; i<m; i++) {
27b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr);
28b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
29b0a32e0cSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr);
300752156aSBarry Smith       }
31b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
32b97cf49bSBarry Smith     }
33b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
34b97cf49bSBarry Smith   }
35b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
363a40ed3dSBarry Smith   PetscFunctionReturn(0);
37b97cf49bSBarry Smith }
38b97cf49bSBarry Smith 
394a2ae208SSatish Balay #undef __FUNCT__
404a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj"
41dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer)
42b97cf49bSBarry Smith {
43dfbe8321SBarry Smith   PetscErrorCode ierr;
4432077d6dSBarry Smith   PetscTruth iascii;
45b97cf49bSBarry Smith 
46433994e6SBarry Smith   PetscFunctionBegin;
4732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
4832077d6dSBarry Smith   if (iascii) {
493eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
505cd90555SBarry Smith   } else {
5129bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name);
52b97cf49bSBarry Smith   }
533a40ed3dSBarry Smith   PetscFunctionReturn(0);
54b97cf49bSBarry Smith }
55b97cf49bSBarry Smith 
564a2ae208SSatish Balay #undef __FUNCT__
574a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj"
58dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat)
59b97cf49bSBarry Smith {
603eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data;
61dfbe8321SBarry Smith   PetscErrorCode ierr;
6294d884c6SBarry Smith 
6394d884c6SBarry Smith   PetscFunctionBegin;
64aa482453SBarry Smith #if defined(PETSC_USE_LOG)
65b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz);
66b97cf49bSBarry Smith #endif
67606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
680400557aSBarry Smith   if (a->freeaij) {
69606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
70606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
711198db28SBarry Smith     if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);}
720400557aSBarry Smith   }
730400557aSBarry Smith   ierr = PetscFree(a->rowners);CHKERRQ(ierr);
741198db28SBarry Smith   ierr = PetscFree(a);CHKERRQ(ierr);
753a40ed3dSBarry Smith   PetscFunctionReturn(0);
76b97cf49bSBarry Smith }
77b97cf49bSBarry Smith 
784a2ae208SSatish Balay #undef __FUNCT__
794a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
80dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op)
81b97cf49bSBarry Smith {
823eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
83b97cf49bSBarry Smith 
84433994e6SBarry Smith   PetscFunctionBegin;
8512c028f9SKris Buschelman   switch (op) {
8677e54ba9SKris Buschelman   case MAT_SYMMETRIC:
8712c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
889a4540c5SBarry Smith   case MAT_HERMITIAN:
89b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
9012c028f9SKris Buschelman     break;
919a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
929a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
939a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
949a4540c5SBarry Smith     a->symmetric = PETSC_FALSE;
959a4540c5SBarry Smith     break;
969a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
979a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
989a4540c5SBarry Smith     break;
9912c028f9SKris Buschelman   default:
100b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
10112c028f9SKris Buschelman     break;
102b97cf49bSBarry Smith   }
1033a40ed3dSBarry Smith   PetscFunctionReturn(0);
104b97cf49bSBarry Smith }
105b97cf49bSBarry Smith 
106b97cf49bSBarry Smith 
107b97cf49bSBarry Smith /*
108b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
109b97cf49bSBarry Smith */
110b97cf49bSBarry Smith 
1114a2ae208SSatish Balay #undef __FUNCT__
1124a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
113dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
114b97cf49bSBarry Smith {
1153eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
116*6849ba73SBarry Smith   PetscErrorCode ierr;
117*6849ba73SBarry Smith   int        i,j,*diag,m = A->m;
118b97cf49bSBarry Smith 
119433994e6SBarry Smith   PetscFunctionBegin;
120b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
121b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
122273d9f13SBarry Smith   for (i=0; i<A->m; i++) {
123b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
124b97cf49bSBarry Smith       if (a->j[j] == i) {
125b97cf49bSBarry Smith         diag[i] = j;
126b97cf49bSBarry Smith         break;
127b97cf49bSBarry Smith       }
128b97cf49bSBarry Smith     }
129b97cf49bSBarry Smith   }
130b97cf49bSBarry Smith   a->diag = diag;
1313a40ed3dSBarry Smith   PetscFunctionReturn(0);
132b97cf49bSBarry Smith }
133b97cf49bSBarry Smith 
1344a2ae208SSatish Balay #undef __FUNCT__
1354a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
136dfbe8321SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v)
137b97cf49bSBarry Smith {
1383eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
139b97cf49bSBarry Smith   int        *itmp;
140b97cf49bSBarry Smith 
141433994e6SBarry Smith   PetscFunctionBegin;
1420752156aSBarry Smith   row -= a->rstart;
1430752156aSBarry Smith 
144273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
145b97cf49bSBarry Smith 
146b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
147b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
148b97cf49bSBarry Smith   if (idx) {
149b97cf49bSBarry Smith     itmp = a->j + a->i[row];
150b97cf49bSBarry Smith     if (*nz) {
151b97cf49bSBarry Smith       *idx = itmp;
152b97cf49bSBarry Smith     }
153b97cf49bSBarry Smith     else *idx = 0;
154b97cf49bSBarry Smith   }
1553a40ed3dSBarry Smith   PetscFunctionReturn(0);
156b97cf49bSBarry Smith }
157b97cf49bSBarry Smith 
1584a2ae208SSatish Balay #undef __FUNCT__
1594a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
160dfbe8321SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v)
161b97cf49bSBarry Smith {
162433994e6SBarry Smith   PetscFunctionBegin;
1633a40ed3dSBarry Smith   PetscFunctionReturn(0);
164b97cf49bSBarry Smith }
165b97cf49bSBarry Smith 
1664a2ae208SSatish Balay #undef __FUNCT__
1674a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj"
168dfbe8321SBarry Smith PetscErrorCode MatGetBlockSize_MPIAdj(Mat A,int *bs)
169b97cf49bSBarry Smith {
170433994e6SBarry Smith   PetscFunctionBegin;
171b97cf49bSBarry Smith   *bs = 1;
1723a40ed3dSBarry Smith   PetscFunctionReturn(0);
173b97cf49bSBarry Smith }
174b97cf49bSBarry Smith 
175b97cf49bSBarry Smith 
1764a2ae208SSatish Balay #undef __FUNCT__
1774a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
178dfbe8321SBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
179b97cf49bSBarry Smith {
1803eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
181dfbe8321SBarry Smith   PetscErrorCode ierr;
1820f5bd95cSBarry Smith   PetscTruth  flag;
183b97cf49bSBarry Smith 
184433994e6SBarry Smith   PetscFunctionBegin;
185b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
186273d9f13SBarry Smith   if ((A->m != B->m) ||(a->nz != b->nz)) {
1870f5bd95cSBarry Smith     flag = PETSC_FALSE;
188b97cf49bSBarry Smith   }
189b97cf49bSBarry Smith 
190b97cf49bSBarry Smith   /* if the a->i are the same */
191273d9f13SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr);
192b97cf49bSBarry Smith 
193b97cf49bSBarry Smith   /* if a->j are the same */
1940f5bd95cSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr);
195b97cf49bSBarry Smith 
196ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
1973a40ed3dSBarry Smith   PetscFunctionReturn(0);
198b97cf49bSBarry Smith }
199b97cf49bSBarry Smith 
2004a2ae208SSatish Balay #undef __FUNCT__
2014a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
202dfbe8321SBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int *ia[],int *ja[],PetscTruth *done)
2036cd91f64SBarry Smith {
204dfbe8321SBarry Smith   PetscErrorCode ierr;
205dfbe8321SBarry Smith   int size,i;
2066cd91f64SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
2076cd91f64SBarry Smith 
2086cd91f64SBarry Smith   PetscFunctionBegin;
2096cd91f64SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2106cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
2116cd91f64SBarry Smith   *m    = A->m;
2126cd91f64SBarry Smith   *ia   = a->i;
2136cd91f64SBarry Smith   *ja   = a->j;
2146cd91f64SBarry Smith   *done = PETSC_TRUE;
215d42735a1SBarry Smith   if (oshift) {
216d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
217d42735a1SBarry Smith       (*ja)[i]++;
218d42735a1SBarry Smith     }
219d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
220d42735a1SBarry Smith   }
221d42735a1SBarry Smith   PetscFunctionReturn(0);
222d42735a1SBarry Smith }
223d42735a1SBarry Smith 
2244a2ae208SSatish Balay #undef __FUNCT__
2254a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
226dfbe8321SBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int *ia[],int *ja[],PetscTruth *done)
227d42735a1SBarry Smith {
228d42735a1SBarry Smith   int        i;
229d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
230d42735a1SBarry Smith 
231d42735a1SBarry Smith   PetscFunctionBegin;
232c5b3d447SBarry Smith   if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()");
233c5b3d447SBarry Smith   if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()");
234d42735a1SBarry Smith   if (oshift) {
235d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
236d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
237d42735a1SBarry Smith       (*ja)[i]--;
238d42735a1SBarry Smith     }
239d42735a1SBarry Smith   }
2406cd91f64SBarry Smith   PetscFunctionReturn(0);
2416cd91f64SBarry Smith }
242b97cf49bSBarry Smith 
243b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2440b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2453eda8832SBarry Smith        MatGetRow_MPIAdj,
2463eda8832SBarry Smith        MatRestoreRow_MPIAdj,
247b97cf49bSBarry Smith        0,
24897304618SKris Buschelman /* 4*/ 0,
249b97cf49bSBarry Smith        0,
250b97cf49bSBarry Smith        0,
251b97cf49bSBarry Smith        0,
2520b3b1487SBarry Smith        0,
2530b3b1487SBarry Smith        0,
25497304618SKris Buschelman /*10*/ 0,
2550b3b1487SBarry Smith        0,
2560b3b1487SBarry Smith        0,
2570b3b1487SBarry Smith        0,
2580b3b1487SBarry Smith        0,
25997304618SKris Buschelman /*15*/ 0,
2603eda8832SBarry Smith        MatEqual_MPIAdj,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        0,
2630b3b1487SBarry Smith        0,
26497304618SKris Buschelman /*20*/ 0,
2650b3b1487SBarry Smith        0,
2660b3b1487SBarry Smith        0,
2673eda8832SBarry Smith        MatSetOption_MPIAdj,
2680b3b1487SBarry Smith        0,
26997304618SKris Buschelman /*25*/ 0,
2700b3b1487SBarry Smith        0,
2710b3b1487SBarry Smith        0,
2720b3b1487SBarry Smith        0,
2730b3b1487SBarry Smith        0,
27497304618SKris Buschelman /*30*/ 0,
2750b3b1487SBarry Smith        0,
276273d9f13SBarry Smith        0,
277273d9f13SBarry Smith        0,
2780b3b1487SBarry Smith        0,
27997304618SKris Buschelman /*35*/ 0,
2800b3b1487SBarry Smith        0,
2810b3b1487SBarry Smith        0,
2820b3b1487SBarry Smith        0,
2830b3b1487SBarry Smith        0,
28497304618SKris Buschelman /*40*/ 0,
2850b3b1487SBarry Smith        0,
2860b3b1487SBarry Smith        0,
2870b3b1487SBarry Smith        0,
2880b3b1487SBarry Smith        0,
28997304618SKris Buschelman /*45*/ 0,
2900b3b1487SBarry Smith        0,
2910b3b1487SBarry Smith        0,
2920b3b1487SBarry Smith        0,
2930b3b1487SBarry Smith        0,
29497304618SKris Buschelman /*50*/ MatGetBlockSize_MPIAdj,
2956cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
296d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
297b97cf49bSBarry Smith        0,
298b97cf49bSBarry Smith        0,
29997304618SKris Buschelman /*55*/ 0,
300b97cf49bSBarry Smith        0,
3010752156aSBarry Smith        0,
3020752156aSBarry Smith        0,
3030b3b1487SBarry Smith        0,
30497304618SKris Buschelman /*60*/ 0,
305b9b97703SBarry Smith        MatDestroy_MPIAdj,
306b9b97703SBarry Smith        MatView_MPIAdj,
30797304618SKris Buschelman        MatGetPetscMaps_Petsc,
30897304618SKris Buschelman        0,
30997304618SKris Buschelman /*65*/ 0,
31097304618SKris Buschelman        0,
31197304618SKris Buschelman        0,
31297304618SKris Buschelman        0,
31397304618SKris Buschelman        0,
31497304618SKris Buschelman /*70*/ 0,
31597304618SKris Buschelman        0,
31697304618SKris Buschelman        0,
31797304618SKris Buschelman        0,
31897304618SKris Buschelman        0,
31997304618SKris Buschelman /*75*/ 0,
32097304618SKris Buschelman        0,
32197304618SKris Buschelman        0,
32297304618SKris Buschelman        0,
32397304618SKris Buschelman        0,
32497304618SKris Buschelman /*80*/ 0,
32597304618SKris Buschelman        0,
32697304618SKris Buschelman        0,
32797304618SKris Buschelman        0,
32897304618SKris Buschelman /*85*/ 0
32997304618SKris Buschelman };
330b97cf49bSBarry Smith 
331a23d5eceSKris Buschelman EXTERN_C_BEGIN
332a23d5eceSKris Buschelman #undef __FUNCT__
333a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
334dfbe8321SBarry Smith PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,int *i,int *j,int *values)
335a23d5eceSKris Buschelman {
336a23d5eceSKris Buschelman   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
337dfbe8321SBarry Smith   PetscErrorCode ierr;
338a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g)
339a23d5eceSKris Buschelman   int        ii;
340a23d5eceSKris Buschelman #endif
341a23d5eceSKris Buschelman 
342a23d5eceSKris Buschelman   PetscFunctionBegin;
343a23d5eceSKris Buschelman   B->preallocated = PETSC_TRUE;
344a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g)
345a23d5eceSKris Buschelman   if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]);
346a23d5eceSKris Buschelman   for (ii=1; ii<B->m; ii++) {
347a23d5eceSKris Buschelman     if (i[ii] < 0 || i[ii] < i[ii-1]) {
348a23d5eceSKris Buschelman       SETERRQ4(1,"i[%d]=%d index is out of range: i[%d]=%d",ii,i[ii],ii-1,i[ii-1]);
349a23d5eceSKris Buschelman     }
350a23d5eceSKris Buschelman   }
351a23d5eceSKris Buschelman   for (ii=0; ii<i[B->m]; ii++) {
352a23d5eceSKris Buschelman     if (j[ii] < 0 || j[ii] >= B->N) {
353a23d5eceSKris Buschelman       SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]);
354a23d5eceSKris Buschelman     }
355a23d5eceSKris Buschelman   }
356a23d5eceSKris Buschelman #endif
357a23d5eceSKris Buschelman 
358a23d5eceSKris Buschelman   b->j      = j;
359a23d5eceSKris Buschelman   b->i      = i;
360a23d5eceSKris Buschelman   b->values = values;
361a23d5eceSKris Buschelman 
362a23d5eceSKris Buschelman   b->nz               = i[B->m];
363a23d5eceSKris Buschelman   b->diag             = 0;
364a23d5eceSKris Buschelman   b->symmetric        = PETSC_FALSE;
365a23d5eceSKris Buschelman   b->freeaij          = PETSC_TRUE;
366a23d5eceSKris Buschelman 
367a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
368a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
369a23d5eceSKris Buschelman   PetscFunctionReturn(0);
370a23d5eceSKris Buschelman }
371a23d5eceSKris Buschelman EXTERN_C_END
372b97cf49bSBarry Smith 
3730bad9183SKris Buschelman /*MC
374fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
3750bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
3760bad9183SKris Buschelman 
3770bad9183SKris Buschelman   Level: beginner
3780bad9183SKris Buschelman 
3790bad9183SKris Buschelman .seealso: MatCreateMPIAdj
3800bad9183SKris Buschelman M*/
3810bad9183SKris Buschelman 
382273d9f13SBarry Smith EXTERN_C_BEGIN
3834a2ae208SSatish Balay #undef __FUNCT__
3844a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
385dfbe8321SBarry Smith PetscErrorCode MatCreate_MPIAdj(Mat B)
386273d9f13SBarry Smith {
387273d9f13SBarry Smith   Mat_MPIAdj *b;
388*6849ba73SBarry Smith   PetscErrorCode ierr;
389*6849ba73SBarry Smith   int        ii,size,rank;
390273d9f13SBarry Smith 
391273d9f13SBarry Smith   PetscFunctionBegin;
392273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
393273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
394273d9f13SBarry Smith 
395b0a32e0cSBarry Smith   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
396b0a32e0cSBarry Smith   B->data             = (void*)b;
397273d9f13SBarry Smith   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
398273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
399273d9f13SBarry Smith   B->factor           = 0;
400273d9f13SBarry Smith   B->lupivotthreshold = 1.0;
401273d9f13SBarry Smith   B->mapping          = 0;
402273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
403273d9f13SBarry Smith 
404273d9f13SBarry Smith   ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr);
405273d9f13SBarry Smith   B->n = B->N;
406273d9f13SBarry Smith 
407273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
408273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
4098a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
410273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
4118a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
412273d9f13SBarry Smith 
413b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr);
414b0a32e0cSBarry Smith   PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
415273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
416273d9f13SBarry Smith   b->rowners[0] = 0;
417273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
418273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
419273d9f13SBarry Smith   }
420273d9f13SBarry Smith   b->rstart = b->rowners[rank];
421273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
422a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C",
423a23d5eceSKris Buschelman                                     "MatMPIAdjSetPreallocation_MPIAdj",
424a23d5eceSKris Buschelman                                      MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
425273d9f13SBarry Smith   PetscFunctionReturn(0);
426273d9f13SBarry Smith }
427273d9f13SBarry Smith EXTERN_C_END
428273d9f13SBarry Smith 
4294a2ae208SSatish Balay #undef __FUNCT__
4304a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
431a23d5eceSKris Buschelman /*@C
432a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
433a23d5eceSKris Buschelman 
434a23d5eceSKris Buschelman    Collective on MPI_Comm
435a23d5eceSKris Buschelman 
436a23d5eceSKris Buschelman    Input Parameters:
437a23d5eceSKris Buschelman +  A - the matrix
438a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
439a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
440a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
441a23d5eceSKris Buschelman -  values - [optional] edge weights
442a23d5eceSKris Buschelman 
443a23d5eceSKris Buschelman    Level: intermediate
444a23d5eceSKris Buschelman 
445a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
446a23d5eceSKris Buschelman @*/
447dfbe8321SBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values)
448273d9f13SBarry Smith {
449dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,int*,int*,int*);
450273d9f13SBarry Smith 
451273d9f13SBarry Smith   PetscFunctionBegin;
452a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
453a23d5eceSKris Buschelman   if (f) {
454a23d5eceSKris Buschelman     ierr = (*f)(B,i,j,values);CHKERRQ(ierr);
455273d9f13SBarry Smith   }
456273d9f13SBarry Smith   PetscFunctionReturn(0);
457273d9f13SBarry Smith }
458273d9f13SBarry Smith 
4594a2ae208SSatish Balay #undef __FUNCT__
4604a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
461b97cf49bSBarry Smith /*@C
4623eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
463b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
464b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
465b97cf49bSBarry Smith 
466ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
467ef5ee4d1SLois Curfman McInnes 
468b97cf49bSBarry Smith    Input Parameters:
469c2e958feSBarry Smith +  comm - MPI communicator
4700752156aSBarry Smith .  m - number of local rows
471b97cf49bSBarry Smith .  n - number of columns
472b97cf49bSBarry Smith .  i - the indices into j for the start of each row
473329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
474ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
475329f5518SBarry Smith -  values -[optional] edge weights
476b97cf49bSBarry Smith 
477b97cf49bSBarry Smith    Output Parameter:
478b97cf49bSBarry Smith .  A - the matrix
479b97cf49bSBarry Smith 
48015091d37SBarry Smith    Level: intermediate
48115091d37SBarry Smith 
4824bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
4834bc6d8c0SBarry Smith    MatSetValues().
484c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
485ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
4861198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
487327e1a73SBarry Smith    Should not include the matrix diagonals.
488b97cf49bSBarry Smith 
489e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
490ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
491ededeb1aSvictorle 
492ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
493b97cf49bSBarry Smith 
494e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
495b97cf49bSBarry Smith @*/
496dfbe8321SBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A)
497b97cf49bSBarry Smith {
498dfbe8321SBarry Smith   PetscErrorCode ierr;
499b97cf49bSBarry Smith 
500433994e6SBarry Smith   PetscFunctionBegin;
501273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
502273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
503273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
5043a40ed3dSBarry Smith   PetscFunctionReturn(0);
505b97cf49bSBarry Smith }
506b97cf49bSBarry Smith 
507273d9f13SBarry Smith EXTERN_C_BEGIN
5084a2ae208SSatish Balay #undef __FUNCT__
5094a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj"
510dfbe8321SBarry Smith PetscErrorCode MatConvertTo_MPIAdj(Mat A,MatType type,Mat *newmat)
511b3cc6726SBarry Smith {
512676c34cdSKris Buschelman   Mat               B;
513*6849ba73SBarry Smith   PetscErrorCode ierr;
514*6849ba73SBarry Smith   int               i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
515b3cc6726SBarry Smith   const int         *rj;
516b3cc6726SBarry Smith   const PetscScalar *ra;
517c2e958feSBarry Smith   MPI_Comm          comm;
518c2e958feSBarry Smith 
519c2e958feSBarry Smith   PetscFunctionBegin;
520c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
521c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
522c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
523c2e958feSBarry Smith 
524c2e958feSBarry Smith   /* count the number of nonzeros per row */
525c2e958feSBarry Smith   for (i=0; i<m; i++) {
526c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
527c2e958feSBarry Smith     for (j=0; j<len; j++) {
528c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
529c2e958feSBarry Smith     }
530c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
531c2e958feSBarry Smith     nzeros += len;
532c2e958feSBarry Smith   }
533c2e958feSBarry Smith 
534c2e958feSBarry Smith   /* malloc space for nonzeros */
535b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr);
536b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr);
537b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr);
538c2e958feSBarry Smith 
539c2e958feSBarry Smith   nzeros = 0;
540c2e958feSBarry Smith   ia[0]  = 0;
541c2e958feSBarry Smith   for (i=0; i<m; i++) {
542c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
543c2e958feSBarry Smith     cnt     = 0;
544c2e958feSBarry Smith     for (j=0; j<len; j++) {
545c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
5464c49b128SBarry Smith         a[nzeros+cnt]    = (int) PetscAbsScalar(ra[j]);
547c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
548c2e958feSBarry Smith       }
549c2e958feSBarry Smith     }
550c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
551c2e958feSBarry Smith     nzeros += cnt;
552c2e958feSBarry Smith     ia[i+1] = nzeros;
553c2e958feSBarry Smith   }
554c2e958feSBarry Smith 
555c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
556f204ca49SKris Buschelman   ierr = MatCreate(comm,m,N,PETSC_DETERMINE,N,&B);CHKERRQ(ierr);
557f204ca49SKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
558f204ca49SKris Buschelman   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
559c2e958feSBarry Smith 
560676c34cdSKris Buschelman   /* Fake support for "inplace" convert. */
561676c34cdSKris Buschelman   if (*newmat == A) {
562676c34cdSKris Buschelman     ierr = MatDestroy(A);CHKERRQ(ierr);
563676c34cdSKris Buschelman   }
564676c34cdSKris Buschelman   *newmat = B;
565c2e958feSBarry Smith   PetscFunctionReturn(0);
566c2e958feSBarry Smith }
567273d9f13SBarry Smith EXTERN_C_END
568433994e6SBarry Smith 
569433994e6SBarry Smith 
570433994e6SBarry Smith 
571433994e6SBarry Smith 
572433994e6SBarry Smith 
573