xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision ae77e1216f63a101c7dbd925cdbaa61ef86e6ca7)
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;
13b24ad042SBarry Smith   PetscInt          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++) {
2777431f27SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+a->rstart);CHKERRQ(ierr);
28b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
2977431f27SBarry 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 {
5179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"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)
6577431f27SBarry 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);
75901853e0SKris Buschelman 
76901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
773a40ed3dSBarry Smith   PetscFunctionReturn(0);
78b97cf49bSBarry Smith }
79b97cf49bSBarry Smith 
804a2ae208SSatish Balay #undef __FUNCT__
814a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
82dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op)
83b97cf49bSBarry Smith {
843eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
85b97cf49bSBarry Smith 
86433994e6SBarry Smith   PetscFunctionBegin;
8712c028f9SKris Buschelman   switch (op) {
8877e54ba9SKris Buschelman   case MAT_SYMMETRIC:
8912c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
909a4540c5SBarry Smith   case MAT_HERMITIAN:
91b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
9212c028f9SKris Buschelman     break;
939a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
949a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
959a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
969a4540c5SBarry Smith     a->symmetric = PETSC_FALSE;
979a4540c5SBarry Smith     break;
989a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
999a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
1009a4540c5SBarry Smith     break;
10112c028f9SKris Buschelman   default:
102b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
10312c028f9SKris Buschelman     break;
104b97cf49bSBarry Smith   }
1053a40ed3dSBarry Smith   PetscFunctionReturn(0);
106b97cf49bSBarry Smith }
107b97cf49bSBarry Smith 
108b97cf49bSBarry Smith 
109b97cf49bSBarry Smith /*
110b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
111b97cf49bSBarry Smith */
112b97cf49bSBarry Smith 
1134a2ae208SSatish Balay #undef __FUNCT__
1144a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
115dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
116b97cf49bSBarry Smith {
1173eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
1186849ba73SBarry Smith   PetscErrorCode ierr;
119b24ad042SBarry Smith   PetscInt       i,j,*diag,m = A->m;
120b97cf49bSBarry Smith 
121433994e6SBarry Smith   PetscFunctionBegin;
122b24ad042SBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&diag);CHKERRQ(ierr);
123b24ad042SBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(PetscInt));
124273d9f13SBarry Smith   for (i=0; i<A->m; i++) {
125b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
126b97cf49bSBarry Smith       if (a->j[j] == i) {
127b97cf49bSBarry Smith         diag[i] = j;
128b97cf49bSBarry Smith         break;
129b97cf49bSBarry Smith       }
130b97cf49bSBarry Smith     }
131b97cf49bSBarry Smith   }
132b97cf49bSBarry Smith   a->diag = diag;
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
134b97cf49bSBarry Smith }
135b97cf49bSBarry Smith 
1364a2ae208SSatish Balay #undef __FUNCT__
1374a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
138b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
139b97cf49bSBarry Smith {
1403eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
141b24ad042SBarry Smith   PetscInt   *itmp;
142b97cf49bSBarry Smith 
143433994e6SBarry Smith   PetscFunctionBegin;
1440752156aSBarry Smith   row -= a->rstart;
1450752156aSBarry Smith 
146273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
147b97cf49bSBarry Smith 
148b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
149b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
150b97cf49bSBarry Smith   if (idx) {
151b97cf49bSBarry Smith     itmp = a->j + a->i[row];
152b97cf49bSBarry Smith     if (*nz) {
153b97cf49bSBarry Smith       *idx = itmp;
154b97cf49bSBarry Smith     }
155b97cf49bSBarry Smith     else *idx = 0;
156b97cf49bSBarry Smith   }
1573a40ed3dSBarry Smith   PetscFunctionReturn(0);
158b97cf49bSBarry Smith }
159b97cf49bSBarry Smith 
1604a2ae208SSatish Balay #undef __FUNCT__
1614a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
162b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
163b97cf49bSBarry Smith {
164433994e6SBarry Smith   PetscFunctionBegin;
1653a40ed3dSBarry Smith   PetscFunctionReturn(0);
166b97cf49bSBarry Smith }
167b97cf49bSBarry Smith 
1684a2ae208SSatish Balay #undef __FUNCT__
1694a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj"
170b24ad042SBarry Smith PetscErrorCode MatGetBlockSize_MPIAdj(Mat A,PetscInt *bs)
171b97cf49bSBarry Smith {
172433994e6SBarry Smith   PetscFunctionBegin;
173b97cf49bSBarry Smith   *bs = 1;
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
175b97cf49bSBarry Smith }
176b97cf49bSBarry Smith 
177b97cf49bSBarry Smith 
1784a2ae208SSatish Balay #undef __FUNCT__
1794a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
180dfbe8321SBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
181b97cf49bSBarry Smith {
1823eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
183dfbe8321SBarry Smith   PetscErrorCode ierr;
1840f5bd95cSBarry Smith   PetscTruth     flag;
185b97cf49bSBarry Smith 
186433994e6SBarry Smith   PetscFunctionBegin;
187b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
188273d9f13SBarry Smith   if ((A->m != B->m) ||(a->nz != b->nz)) {
1890f5bd95cSBarry Smith     flag = PETSC_FALSE;
190b97cf49bSBarry Smith   }
191b97cf49bSBarry Smith 
192b97cf49bSBarry Smith   /* if the a->i are the same */
193b24ad042SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
194b97cf49bSBarry Smith 
195b97cf49bSBarry Smith   /* if a->j are the same */
196b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
197b97cf49bSBarry Smith 
198ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
1993a40ed3dSBarry Smith   PetscFunctionReturn(0);
200b97cf49bSBarry Smith }
201b97cf49bSBarry Smith 
2024a2ae208SSatish Balay #undef __FUNCT__
2034a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
204b24ad042SBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
2056cd91f64SBarry Smith {
206dfbe8321SBarry Smith   PetscErrorCode ierr;
207b24ad042SBarry Smith   PetscMPIInt    size;
208b24ad042SBarry Smith   PetscInt       i;
2096cd91f64SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data;
2106cd91f64SBarry Smith 
2116cd91f64SBarry Smith   PetscFunctionBegin;
2126cd91f64SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2136cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
2146cd91f64SBarry Smith   *m    = A->m;
2156cd91f64SBarry Smith   *ia   = a->i;
2166cd91f64SBarry Smith   *ja   = a->j;
2176cd91f64SBarry Smith   *done = PETSC_TRUE;
218d42735a1SBarry Smith   if (oshift) {
219d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
220d42735a1SBarry Smith       (*ja)[i]++;
221d42735a1SBarry Smith     }
222d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
223d42735a1SBarry Smith   }
224d42735a1SBarry Smith   PetscFunctionReturn(0);
225d42735a1SBarry Smith }
226d42735a1SBarry Smith 
2274a2ae208SSatish Balay #undef __FUNCT__
2284a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
229b24ad042SBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
230d42735a1SBarry Smith {
231b24ad042SBarry Smith   PetscInt   i;
232d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
233d42735a1SBarry Smith 
234d42735a1SBarry Smith   PetscFunctionBegin;
235e005ede5SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
236e005ede5SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
237d42735a1SBarry Smith   if (oshift) {
238d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
239d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
240d42735a1SBarry Smith       (*ja)[i]--;
241d42735a1SBarry Smith     }
242d42735a1SBarry Smith   }
2436cd91f64SBarry Smith   PetscFunctionReturn(0);
2446cd91f64SBarry Smith }
245b97cf49bSBarry Smith 
246b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2470b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2483eda8832SBarry Smith        MatGetRow_MPIAdj,
2493eda8832SBarry Smith        MatRestoreRow_MPIAdj,
250b97cf49bSBarry Smith        0,
25197304618SKris Buschelman /* 4*/ 0,
252b97cf49bSBarry Smith        0,
253b97cf49bSBarry Smith        0,
254b97cf49bSBarry Smith        0,
2550b3b1487SBarry Smith        0,
2560b3b1487SBarry Smith        0,
25797304618SKris Buschelman /*10*/ 0,
2580b3b1487SBarry Smith        0,
2590b3b1487SBarry Smith        0,
2600b3b1487SBarry Smith        0,
2610b3b1487SBarry Smith        0,
26297304618SKris Buschelman /*15*/ 0,
2633eda8832SBarry Smith        MatEqual_MPIAdj,
2640b3b1487SBarry Smith        0,
2650b3b1487SBarry Smith        0,
2660b3b1487SBarry Smith        0,
26797304618SKris Buschelman /*20*/ 0,
2680b3b1487SBarry Smith        0,
2690b3b1487SBarry Smith        0,
2703eda8832SBarry Smith        MatSetOption_MPIAdj,
2710b3b1487SBarry Smith        0,
27297304618SKris Buschelman /*25*/ 0,
2730b3b1487SBarry Smith        0,
2740b3b1487SBarry Smith        0,
2750b3b1487SBarry Smith        0,
2760b3b1487SBarry Smith        0,
27797304618SKris Buschelman /*30*/ 0,
2780b3b1487SBarry Smith        0,
279273d9f13SBarry Smith        0,
280273d9f13SBarry Smith        0,
2810b3b1487SBarry Smith        0,
28297304618SKris Buschelman /*35*/ 0,
2830b3b1487SBarry Smith        0,
2840b3b1487SBarry Smith        0,
2850b3b1487SBarry Smith        0,
2860b3b1487SBarry Smith        0,
28797304618SKris Buschelman /*40*/ 0,
2880b3b1487SBarry Smith        0,
2890b3b1487SBarry Smith        0,
2900b3b1487SBarry Smith        0,
2910b3b1487SBarry Smith        0,
29297304618SKris Buschelman /*45*/ 0,
2930b3b1487SBarry Smith        0,
2940b3b1487SBarry Smith        0,
2950b3b1487SBarry Smith        0,
2960b3b1487SBarry Smith        0,
29797304618SKris Buschelman /*50*/ MatGetBlockSize_MPIAdj,
2986cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
299d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
300b97cf49bSBarry Smith        0,
301b97cf49bSBarry Smith        0,
30297304618SKris Buschelman /*55*/ 0,
303b97cf49bSBarry Smith        0,
3040752156aSBarry Smith        0,
3050752156aSBarry Smith        0,
3060b3b1487SBarry Smith        0,
30797304618SKris Buschelman /*60*/ 0,
308b9b97703SBarry Smith        MatDestroy_MPIAdj,
309b9b97703SBarry Smith        MatView_MPIAdj,
31097304618SKris Buschelman        MatGetPetscMaps_Petsc,
31197304618SKris Buschelman        0,
31297304618SKris Buschelman /*65*/ 0,
31397304618SKris Buschelman        0,
31497304618SKris Buschelman        0,
31597304618SKris Buschelman        0,
31697304618SKris Buschelman        0,
31797304618SKris Buschelman /*70*/ 0,
31897304618SKris Buschelman        0,
31997304618SKris Buschelman        0,
32097304618SKris Buschelman        0,
32197304618SKris Buschelman        0,
32297304618SKris Buschelman /*75*/ 0,
32397304618SKris Buschelman        0,
32497304618SKris Buschelman        0,
32597304618SKris Buschelman        0,
32697304618SKris Buschelman        0,
32797304618SKris Buschelman /*80*/ 0,
32897304618SKris Buschelman        0,
32997304618SKris Buschelman        0,
33097304618SKris Buschelman        0,
331865e5f61SKris Buschelman        0,
332865e5f61SKris Buschelman /*85*/ 0,
333865e5f61SKris Buschelman        0,
334865e5f61SKris Buschelman        0,
335865e5f61SKris Buschelman        0,
336865e5f61SKris Buschelman        0,
337865e5f61SKris Buschelman /*90*/ 0,
338865e5f61SKris Buschelman        0,
339865e5f61SKris Buschelman        0,
340865e5f61SKris Buschelman        0,
341865e5f61SKris Buschelman        0,
342865e5f61SKris Buschelman /*95*/ 0,
343865e5f61SKris Buschelman        0,
344865e5f61SKris Buschelman        0,
345865e5f61SKris Buschelman        0};
346b97cf49bSBarry Smith 
347a23d5eceSKris Buschelman EXTERN_C_BEGIN
348a23d5eceSKris Buschelman #undef __FUNCT__
349a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
350b24ad042SBarry Smith PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
351a23d5eceSKris Buschelman {
352a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj *)B->data;
353dfbe8321SBarry Smith   PetscErrorCode ierr;
354a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g)
355b24ad042SBarry Smith   PetscInt       ii;
356a23d5eceSKris Buschelman #endif
357a23d5eceSKris Buschelman 
358a23d5eceSKris Buschelman   PetscFunctionBegin;
359a23d5eceSKris Buschelman   B->preallocated = PETSC_TRUE;
360a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g)
36177431f27SBarry Smith   if (i[0] != 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]);
362a23d5eceSKris Buschelman   for (ii=1; ii<B->m; ii++) {
363a23d5eceSKris Buschelman     if (i[ii] < 0 || i[ii] < i[ii-1]) {
36477431f27SBarry Smith       SETERRQ4(PETSC_ERR_ARG_OUTOFRANGE,"i[%D]=%D index is out of range: i[%D]=%D",ii,i[ii],ii-1,i[ii-1]);
365a23d5eceSKris Buschelman     }
366a23d5eceSKris Buschelman   }
367a23d5eceSKris Buschelman   for (ii=0; ii<i[B->m]; ii++) {
368a23d5eceSKris Buschelman     if (j[ii] < 0 || j[ii] >= B->N) {
36977431f27SBarry Smith       SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]);
370a23d5eceSKris Buschelman     }
371a23d5eceSKris Buschelman   }
372a23d5eceSKris Buschelman #endif
373a23d5eceSKris Buschelman 
374a23d5eceSKris Buschelman   b->j      = j;
375a23d5eceSKris Buschelman   b->i      = i;
376a23d5eceSKris Buschelman   b->values = values;
377a23d5eceSKris Buschelman 
378a23d5eceSKris Buschelman   b->nz               = i[B->m];
379a23d5eceSKris Buschelman   b->diag             = 0;
380a23d5eceSKris Buschelman   b->symmetric        = PETSC_FALSE;
381a23d5eceSKris Buschelman   b->freeaij          = PETSC_TRUE;
382a23d5eceSKris Buschelman 
383a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
384a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
385a23d5eceSKris Buschelman   PetscFunctionReturn(0);
386a23d5eceSKris Buschelman }
387a23d5eceSKris Buschelman EXTERN_C_END
388b97cf49bSBarry Smith 
3890bad9183SKris Buschelman /*MC
390fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
3910bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
3920bad9183SKris Buschelman 
3930bad9183SKris Buschelman   Level: beginner
3940bad9183SKris Buschelman 
3950bad9183SKris Buschelman .seealso: MatCreateMPIAdj
3960bad9183SKris Buschelman M*/
3970bad9183SKris Buschelman 
398273d9f13SBarry Smith EXTERN_C_BEGIN
3994a2ae208SSatish Balay #undef __FUNCT__
4004a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
401dfbe8321SBarry Smith PetscErrorCode MatCreate_MPIAdj(Mat B)
402273d9f13SBarry Smith {
403273d9f13SBarry Smith   Mat_MPIAdj     *b;
4046849ba73SBarry Smith   PetscErrorCode ierr;
405b24ad042SBarry Smith   PetscInt       ii;
406b24ad042SBarry Smith   PetscMPIInt    size,rank;
407273d9f13SBarry Smith 
408273d9f13SBarry Smith   PetscFunctionBegin;
409273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
410273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
411273d9f13SBarry Smith 
412b0a32e0cSBarry Smith   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
413b0a32e0cSBarry Smith   B->data             = (void*)b;
414273d9f13SBarry Smith   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
415273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
416273d9f13SBarry Smith   B->factor           = 0;
417273d9f13SBarry Smith   B->lupivotthreshold = 1.0;
418273d9f13SBarry Smith   B->mapping          = 0;
419273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
420273d9f13SBarry Smith 
421*ae77e121SBarry Smith   ierr = PetscSplitOwnership(B->comm,&B->m,&B->M);CHKERRQ(ierr);
422*ae77e121SBarry Smith   B->n = B->N = PetscMax(B->N,B->n);CHKERRQ(ierr);
423273d9f13SBarry Smith 
424273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
425273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
4268a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
427273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
4288a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
429273d9f13SBarry Smith 
430b24ad042SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&b->rowners);CHKERRQ(ierr);
431b24ad042SBarry Smith   PetscLogObjectMemory(B,(size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
432273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
433273d9f13SBarry Smith   b->rowners[0] = 0;
434273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
435273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
436273d9f13SBarry Smith   }
437273d9f13SBarry Smith   b->rstart = b->rowners[rank];
438273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
439a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C",
440a23d5eceSKris Buschelman                                     "MatMPIAdjSetPreallocation_MPIAdj",
441a23d5eceSKris Buschelman                                      MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
442273d9f13SBarry Smith   PetscFunctionReturn(0);
443273d9f13SBarry Smith }
444273d9f13SBarry Smith EXTERN_C_END
445273d9f13SBarry Smith 
4464a2ae208SSatish Balay #undef __FUNCT__
4474a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
448a23d5eceSKris Buschelman /*@C
449a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
450a23d5eceSKris Buschelman 
451a23d5eceSKris Buschelman    Collective on MPI_Comm
452a23d5eceSKris Buschelman 
453a23d5eceSKris Buschelman    Input Parameters:
454a23d5eceSKris Buschelman +  A - the matrix
455a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
456a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
457a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
458a23d5eceSKris Buschelman -  values - [optional] edge weights
459a23d5eceSKris Buschelman 
460a23d5eceSKris Buschelman    Level: intermediate
461a23d5eceSKris Buschelman 
462a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
463a23d5eceSKris Buschelman @*/
464b24ad042SBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
465273d9f13SBarry Smith {
466b24ad042SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*,PetscInt*);
467273d9f13SBarry Smith 
468273d9f13SBarry Smith   PetscFunctionBegin;
469a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
470a23d5eceSKris Buschelman   if (f) {
471a23d5eceSKris Buschelman     ierr = (*f)(B,i,j,values);CHKERRQ(ierr);
472273d9f13SBarry Smith   }
473273d9f13SBarry Smith   PetscFunctionReturn(0);
474273d9f13SBarry Smith }
475273d9f13SBarry Smith 
4764a2ae208SSatish Balay #undef __FUNCT__
4774a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
478b97cf49bSBarry Smith /*@C
4793eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
480b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
481b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
482b97cf49bSBarry Smith 
483ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
484ef5ee4d1SLois Curfman McInnes 
485b97cf49bSBarry Smith    Input Parameters:
486c2e958feSBarry Smith +  comm - MPI communicator
4870752156aSBarry Smith .  m - number of local rows
488b97cf49bSBarry Smith .  n - number of columns
489b97cf49bSBarry Smith .  i - the indices into j for the start of each row
490329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
491ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
492329f5518SBarry Smith -  values -[optional] edge weights
493b97cf49bSBarry Smith 
494b97cf49bSBarry Smith    Output Parameter:
495b97cf49bSBarry Smith .  A - the matrix
496b97cf49bSBarry Smith 
49715091d37SBarry Smith    Level: intermediate
49815091d37SBarry Smith 
4994bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
5004bc6d8c0SBarry Smith    MatSetValues().
501c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
502ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
5031198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
504327e1a73SBarry Smith    Should not include the matrix diagonals.
505b97cf49bSBarry Smith 
506e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
507ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
508ededeb1aSvictorle 
509ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
510b97cf49bSBarry Smith 
511e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
512b97cf49bSBarry Smith @*/
513b24ad042SBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
514b97cf49bSBarry Smith {
515dfbe8321SBarry Smith   PetscErrorCode ierr;
516b97cf49bSBarry Smith 
517433994e6SBarry Smith   PetscFunctionBegin;
518273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
519273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
520273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
5213a40ed3dSBarry Smith   PetscFunctionReturn(0);
522b97cf49bSBarry Smith }
523b97cf49bSBarry Smith 
524273d9f13SBarry Smith EXTERN_C_BEGIN
5254a2ae208SSatish Balay #undef __FUNCT__
5264a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj"
527dfbe8321SBarry Smith PetscErrorCode MatConvertTo_MPIAdj(Mat A,MatType type,Mat *newmat)
528b3cc6726SBarry Smith {
529676c34cdSKris Buschelman   Mat               B;
5306849ba73SBarry Smith   PetscErrorCode    ierr;
531b24ad042SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
532b24ad042SBarry Smith   const PetscInt    *rj;
533b3cc6726SBarry Smith   const PetscScalar *ra;
534c2e958feSBarry Smith   MPI_Comm          comm;
535c2e958feSBarry Smith 
536c2e958feSBarry Smith   PetscFunctionBegin;
537c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
538c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
539c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
540c2e958feSBarry Smith 
541c2e958feSBarry Smith   /* count the number of nonzeros per row */
542c2e958feSBarry Smith   for (i=0; i<m; i++) {
543c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
544c2e958feSBarry Smith     for (j=0; j<len; j++) {
545c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
546c2e958feSBarry Smith     }
547c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
548c2e958feSBarry Smith     nzeros += len;
549c2e958feSBarry Smith   }
550c2e958feSBarry Smith 
551c2e958feSBarry Smith   /* malloc space for nonzeros */
552b24ad042SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr);
553b24ad042SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
554b24ad042SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr);
555c2e958feSBarry Smith 
556c2e958feSBarry Smith   nzeros = 0;
557c2e958feSBarry Smith   ia[0]  = 0;
558c2e958feSBarry Smith   for (i=0; i<m; i++) {
559c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
560c2e958feSBarry Smith     cnt     = 0;
561c2e958feSBarry Smith     for (j=0; j<len; j++) {
562c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
563b24ad042SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
564c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
565c2e958feSBarry Smith       }
566c2e958feSBarry Smith     }
567c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
568c2e958feSBarry Smith     nzeros += cnt;
569c2e958feSBarry Smith     ia[i+1] = nzeros;
570c2e958feSBarry Smith   }
571c2e958feSBarry Smith 
572c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
573f204ca49SKris Buschelman   ierr = MatCreate(comm,m,N,PETSC_DETERMINE,N,&B);CHKERRQ(ierr);
574f204ca49SKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
575f204ca49SKris Buschelman   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
576c2e958feSBarry Smith 
577676c34cdSKris Buschelman   /* Fake support for "inplace" convert. */
578676c34cdSKris Buschelman   if (*newmat == A) {
579676c34cdSKris Buschelman     ierr = MatDestroy(A);CHKERRQ(ierr);
580676c34cdSKris Buschelman   }
581676c34cdSKris Buschelman   *newmat = B;
582c2e958feSBarry Smith   PetscFunctionReturn(0);
583c2e958feSBarry Smith }
584273d9f13SBarry Smith EXTERN_C_END
585433994e6SBarry Smith 
586433994e6SBarry Smith 
587433994e6SBarry Smith 
588433994e6SBarry Smith 
589433994e6SBarry Smith 
590