xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
2*be1d678aSKris Buschelman 
3b97cf49bSBarry Smith /*
4b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
5b97cf49bSBarry Smith */
63eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h"
7325e03aeSBarry Smith #include "petscsys.h"
8b97cf49bSBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII"
11dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
12b97cf49bSBarry Smith {
133eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
14dfbe8321SBarry Smith   PetscErrorCode    ierr;
15b24ad042SBarry Smith   PetscInt          i,j,m = A->m;
16fb9695e5SSatish Balay   char              *name;
17ce1411ecSBarry Smith   PetscViewerFormat format;
18b97cf49bSBarry Smith 
19433994e6SBarry Smith   PetscFunctionBegin;
203a7fca6bSBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
21b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
22fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
233a40ed3dSBarry Smith     PetscFunctionReturn(0);
24fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
25ffac6cdbSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
260752156aSBarry Smith   } else {
27b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
28b97cf49bSBarry Smith     for (i=0; i<m; i++) {
2977431f27SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+a->rstart);CHKERRQ(ierr);
30b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
3177431f27SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr);
320752156aSBarry Smith       }
33b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
34b97cf49bSBarry Smith     }
35b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
36b97cf49bSBarry Smith   }
37b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
383a40ed3dSBarry Smith   PetscFunctionReturn(0);
39b97cf49bSBarry Smith }
40b97cf49bSBarry Smith 
414a2ae208SSatish Balay #undef __FUNCT__
424a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj"
43dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer)
44b97cf49bSBarry Smith {
45dfbe8321SBarry Smith   PetscErrorCode ierr;
4632077d6dSBarry Smith   PetscTruth     iascii;
47b97cf49bSBarry Smith 
48433994e6SBarry Smith   PetscFunctionBegin;
4932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
5032077d6dSBarry Smith   if (iascii) {
513eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
525cd90555SBarry Smith   } else {
5379a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name);
54b97cf49bSBarry Smith   }
553a40ed3dSBarry Smith   PetscFunctionReturn(0);
56b97cf49bSBarry Smith }
57b97cf49bSBarry Smith 
584a2ae208SSatish Balay #undef __FUNCT__
594a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj"
60dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat)
61b97cf49bSBarry Smith {
623eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)mat->data;
63dfbe8321SBarry Smith   PetscErrorCode ierr;
6494d884c6SBarry Smith 
6594d884c6SBarry Smith   PetscFunctionBegin;
66aa482453SBarry Smith #if defined(PETSC_USE_LOG)
6777431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->m,mat->n,a->nz);
68b97cf49bSBarry Smith #endif
69606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
700400557aSBarry Smith   if (a->freeaij) {
71606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
72606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
731198db28SBarry Smith     if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);}
740400557aSBarry Smith   }
750400557aSBarry Smith   ierr = PetscFree(a->rowners);CHKERRQ(ierr);
761198db28SBarry Smith   ierr = PetscFree(a);CHKERRQ(ierr);
77901853e0SKris Buschelman 
78901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
793a40ed3dSBarry Smith   PetscFunctionReturn(0);
80b97cf49bSBarry Smith }
81b97cf49bSBarry Smith 
824a2ae208SSatish Balay #undef __FUNCT__
834a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
84dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op)
85b97cf49bSBarry Smith {
863eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
8763ba0a88SBarry Smith   PetscErrorCode ierr;
88b97cf49bSBarry Smith 
89433994e6SBarry Smith   PetscFunctionBegin;
9012c028f9SKris Buschelman   switch (op) {
9177e54ba9SKris Buschelman   case MAT_SYMMETRIC:
9212c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
939a4540c5SBarry Smith   case MAT_HERMITIAN:
94b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
9512c028f9SKris Buschelman     break;
969a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
979a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
989a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
999a4540c5SBarry Smith     a->symmetric = PETSC_FALSE;
1009a4540c5SBarry Smith     break;
1019a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1029a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
1039a4540c5SBarry Smith     break;
10412c028f9SKris Buschelman   default:
10563ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatSetOption_MPIAdj:Option ignored\n"));CHKERRQ(ierr);
10612c028f9SKris Buschelman     break;
107b97cf49bSBarry Smith   }
1083a40ed3dSBarry Smith   PetscFunctionReturn(0);
109b97cf49bSBarry Smith }
110b97cf49bSBarry Smith 
111b97cf49bSBarry Smith 
112b97cf49bSBarry Smith /*
113b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
114b97cf49bSBarry Smith */
115b97cf49bSBarry Smith 
1164a2ae208SSatish Balay #undef __FUNCT__
1174a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
118dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
119b97cf49bSBarry Smith {
1203eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
1216849ba73SBarry Smith   PetscErrorCode ierr;
122b24ad042SBarry Smith   PetscInt       i,j,*diag,m = A->m;
123b97cf49bSBarry Smith 
124433994e6SBarry Smith   PetscFunctionBegin;
125b24ad042SBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&diag);CHKERRQ(ierr);
12652e6d16bSBarry Smith   ierr = PetscLogObjectMemory(A,(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
127273d9f13SBarry Smith   for (i=0; i<A->m; i++) {
128b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
129b97cf49bSBarry Smith       if (a->j[j] == i) {
130b97cf49bSBarry Smith         diag[i] = j;
131b97cf49bSBarry Smith         break;
132b97cf49bSBarry Smith       }
133b97cf49bSBarry Smith     }
134b97cf49bSBarry Smith   }
135b97cf49bSBarry Smith   a->diag = diag;
1363a40ed3dSBarry Smith   PetscFunctionReturn(0);
137b97cf49bSBarry Smith }
138b97cf49bSBarry Smith 
1394a2ae208SSatish Balay #undef __FUNCT__
1404a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
141b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
142b97cf49bSBarry Smith {
1433eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
144b24ad042SBarry Smith   PetscInt   *itmp;
145b97cf49bSBarry Smith 
146433994e6SBarry Smith   PetscFunctionBegin;
1470752156aSBarry Smith   row -= a->rstart;
1480752156aSBarry Smith 
149273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
150b97cf49bSBarry Smith 
151b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
152b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
153b97cf49bSBarry Smith   if (idx) {
154b97cf49bSBarry Smith     itmp = a->j + a->i[row];
155b97cf49bSBarry Smith     if (*nz) {
156b97cf49bSBarry Smith       *idx = itmp;
157b97cf49bSBarry Smith     }
158b97cf49bSBarry Smith     else *idx = 0;
159b97cf49bSBarry Smith   }
1603a40ed3dSBarry Smith   PetscFunctionReturn(0);
161b97cf49bSBarry Smith }
162b97cf49bSBarry Smith 
1634a2ae208SSatish Balay #undef __FUNCT__
1644a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
165b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
166b97cf49bSBarry Smith {
167433994e6SBarry Smith   PetscFunctionBegin;
1683a40ed3dSBarry Smith   PetscFunctionReturn(0);
169b97cf49bSBarry Smith }
170b97cf49bSBarry Smith 
1714a2ae208SSatish Balay #undef __FUNCT__
1724a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
173dfbe8321SBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
174b97cf49bSBarry Smith {
1753eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
176dfbe8321SBarry Smith   PetscErrorCode ierr;
1770f5bd95cSBarry Smith   PetscTruth     flag;
178b97cf49bSBarry Smith 
179433994e6SBarry Smith   PetscFunctionBegin;
180b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
181273d9f13SBarry Smith   if ((A->m != B->m) ||(a->nz != b->nz)) {
1820f5bd95cSBarry Smith     flag = PETSC_FALSE;
183b97cf49bSBarry Smith   }
184b97cf49bSBarry Smith 
185b97cf49bSBarry Smith   /* if the a->i are the same */
186b24ad042SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
187b97cf49bSBarry Smith 
188b97cf49bSBarry Smith   /* if a->j are the same */
189b24ad042SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr);
190b97cf49bSBarry Smith 
191ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
1923a40ed3dSBarry Smith   PetscFunctionReturn(0);
193b97cf49bSBarry Smith }
194b97cf49bSBarry Smith 
1954a2ae208SSatish Balay #undef __FUNCT__
1964a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
197b24ad042SBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
1986cd91f64SBarry Smith {
199dfbe8321SBarry Smith   PetscErrorCode ierr;
200b24ad042SBarry Smith   PetscMPIInt    size;
201b24ad042SBarry Smith   PetscInt       i;
2026cd91f64SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data;
2036cd91f64SBarry Smith 
2046cd91f64SBarry Smith   PetscFunctionBegin;
2056cd91f64SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2066cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
2076cd91f64SBarry Smith   *m    = A->m;
2086cd91f64SBarry Smith   *ia   = a->i;
2096cd91f64SBarry Smith   *ja   = a->j;
2106cd91f64SBarry Smith   *done = PETSC_TRUE;
211d42735a1SBarry Smith   if (oshift) {
212d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
213d42735a1SBarry Smith       (*ja)[i]++;
214d42735a1SBarry Smith     }
215d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
216d42735a1SBarry Smith   }
217d42735a1SBarry Smith   PetscFunctionReturn(0);
218d42735a1SBarry Smith }
219d42735a1SBarry Smith 
2204a2ae208SSatish Balay #undef __FUNCT__
2214a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
222b24ad042SBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
223d42735a1SBarry Smith {
224b24ad042SBarry Smith   PetscInt   i;
225d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
226d42735a1SBarry Smith 
227d42735a1SBarry Smith   PetscFunctionBegin;
228e005ede5SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
229e005ede5SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
230d42735a1SBarry Smith   if (oshift) {
231d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
232d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
233d42735a1SBarry Smith       (*ja)[i]--;
234d42735a1SBarry Smith     }
235d42735a1SBarry Smith   }
2366cd91f64SBarry Smith   PetscFunctionReturn(0);
2376cd91f64SBarry Smith }
238b97cf49bSBarry Smith 
239b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2400b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2413eda8832SBarry Smith        MatGetRow_MPIAdj,
2423eda8832SBarry Smith        MatRestoreRow_MPIAdj,
243b97cf49bSBarry Smith        0,
24497304618SKris Buschelman /* 4*/ 0,
245b97cf49bSBarry Smith        0,
246b97cf49bSBarry Smith        0,
247b97cf49bSBarry Smith        0,
2480b3b1487SBarry Smith        0,
2490b3b1487SBarry Smith        0,
25097304618SKris Buschelman /*10*/ 0,
2510b3b1487SBarry Smith        0,
2520b3b1487SBarry Smith        0,
2530b3b1487SBarry Smith        0,
2540b3b1487SBarry Smith        0,
25597304618SKris Buschelman /*15*/ 0,
2563eda8832SBarry Smith        MatEqual_MPIAdj,
2570b3b1487SBarry Smith        0,
2580b3b1487SBarry Smith        0,
2590b3b1487SBarry Smith        0,
26097304618SKris Buschelman /*20*/ 0,
2610b3b1487SBarry Smith        0,
2620b3b1487SBarry Smith        0,
2633eda8832SBarry Smith        MatSetOption_MPIAdj,
2640b3b1487SBarry Smith        0,
26597304618SKris Buschelman /*25*/ 0,
2660b3b1487SBarry Smith        0,
2670b3b1487SBarry Smith        0,
2680b3b1487SBarry Smith        0,
2690b3b1487SBarry Smith        0,
27097304618SKris Buschelman /*30*/ 0,
2710b3b1487SBarry Smith        0,
272273d9f13SBarry Smith        0,
273273d9f13SBarry Smith        0,
2740b3b1487SBarry Smith        0,
27597304618SKris Buschelman /*35*/ 0,
2760b3b1487SBarry Smith        0,
2770b3b1487SBarry Smith        0,
2780b3b1487SBarry Smith        0,
2790b3b1487SBarry Smith        0,
28097304618SKris Buschelman /*40*/ 0,
2810b3b1487SBarry Smith        0,
2820b3b1487SBarry Smith        0,
2830b3b1487SBarry Smith        0,
2840b3b1487SBarry Smith        0,
28597304618SKris Buschelman /*45*/ 0,
2860b3b1487SBarry Smith        0,
2870b3b1487SBarry Smith        0,
2880b3b1487SBarry Smith        0,
2890b3b1487SBarry Smith        0,
290521d7252SBarry Smith /*50*/ 0,
2916cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
292d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
293b97cf49bSBarry Smith        0,
294b97cf49bSBarry Smith        0,
29597304618SKris Buschelman /*55*/ 0,
296b97cf49bSBarry Smith        0,
2970752156aSBarry Smith        0,
2980752156aSBarry Smith        0,
2990b3b1487SBarry Smith        0,
30097304618SKris Buschelman /*60*/ 0,
301b9b97703SBarry Smith        MatDestroy_MPIAdj,
302b9b97703SBarry Smith        MatView_MPIAdj,
30397304618SKris Buschelman        MatGetPetscMaps_Petsc,
30497304618SKris Buschelman        0,
30597304618SKris Buschelman /*65*/ 0,
30697304618SKris Buschelman        0,
30797304618SKris Buschelman        0,
30897304618SKris Buschelman        0,
30997304618SKris Buschelman        0,
31097304618SKris Buschelman /*70*/ 0,
31197304618SKris Buschelman        0,
31297304618SKris Buschelman        0,
31397304618SKris Buschelman        0,
31497304618SKris Buschelman        0,
31597304618SKris Buschelman /*75*/ 0,
31697304618SKris Buschelman        0,
31797304618SKris Buschelman        0,
31897304618SKris Buschelman        0,
31997304618SKris Buschelman        0,
32097304618SKris Buschelman /*80*/ 0,
32197304618SKris Buschelman        0,
32297304618SKris Buschelman        0,
32397304618SKris Buschelman        0,
324865e5f61SKris Buschelman        0,
325865e5f61SKris Buschelman /*85*/ 0,
326865e5f61SKris Buschelman        0,
327865e5f61SKris Buschelman        0,
328865e5f61SKris Buschelman        0,
329865e5f61SKris Buschelman        0,
330865e5f61SKris Buschelman /*90*/ 0,
331865e5f61SKris Buschelman        0,
332865e5f61SKris Buschelman        0,
333865e5f61SKris Buschelman        0,
334865e5f61SKris Buschelman        0,
335865e5f61SKris Buschelman /*95*/ 0,
336865e5f61SKris Buschelman        0,
337865e5f61SKris Buschelman        0,
338865e5f61SKris Buschelman        0};
339b97cf49bSBarry Smith 
340a23d5eceSKris Buschelman EXTERN_C_BEGIN
341a23d5eceSKris Buschelman #undef __FUNCT__
342a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
343*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
344a23d5eceSKris Buschelman {
345a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj *)B->data;
346dfbe8321SBarry Smith   PetscErrorCode ierr;
3472515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
348b24ad042SBarry Smith   PetscInt       ii;
349a23d5eceSKris Buschelman #endif
350a23d5eceSKris Buschelman 
351a23d5eceSKris Buschelman   PetscFunctionBegin;
352a23d5eceSKris Buschelman   B->preallocated = PETSC_TRUE;
3532515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
35477431f27SBarry Smith   if (i[0] != 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]);
355a23d5eceSKris Buschelman   for (ii=1; ii<B->m; ii++) {
356a23d5eceSKris Buschelman     if (i[ii] < 0 || i[ii] < i[ii-1]) {
35777431f27SBarry 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]);
358a23d5eceSKris Buschelman     }
359a23d5eceSKris Buschelman   }
360a23d5eceSKris Buschelman   for (ii=0; ii<i[B->m]; ii++) {
361a23d5eceSKris Buschelman     if (j[ii] < 0 || j[ii] >= B->N) {
36277431f27SBarry Smith       SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]);
363a23d5eceSKris Buschelman     }
364a23d5eceSKris Buschelman   }
365a23d5eceSKris Buschelman #endif
366a23d5eceSKris Buschelman 
367a23d5eceSKris Buschelman   b->j      = j;
368a23d5eceSKris Buschelman   b->i      = i;
369a23d5eceSKris Buschelman   b->values = values;
370a23d5eceSKris Buschelman 
371a23d5eceSKris Buschelman   b->nz               = i[B->m];
372a23d5eceSKris Buschelman   b->diag             = 0;
373a23d5eceSKris Buschelman   b->symmetric        = PETSC_FALSE;
374a23d5eceSKris Buschelman   b->freeaij          = PETSC_TRUE;
375a23d5eceSKris Buschelman 
376a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
377a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
378a23d5eceSKris Buschelman   PetscFunctionReturn(0);
379a23d5eceSKris Buschelman }
380a23d5eceSKris Buschelman EXTERN_C_END
381b97cf49bSBarry Smith 
3820bad9183SKris Buschelman /*MC
383fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
3840bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
3850bad9183SKris Buschelman 
3860bad9183SKris Buschelman   Level: beginner
3870bad9183SKris Buschelman 
3880bad9183SKris Buschelman .seealso: MatCreateMPIAdj
3890bad9183SKris Buschelman M*/
3900bad9183SKris Buschelman 
391273d9f13SBarry Smith EXTERN_C_BEGIN
3924a2ae208SSatish Balay #undef __FUNCT__
3934a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
394*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIAdj(Mat B)
395273d9f13SBarry Smith {
396273d9f13SBarry Smith   Mat_MPIAdj     *b;
3976849ba73SBarry Smith   PetscErrorCode ierr;
398b24ad042SBarry Smith   PetscInt       ii;
399b24ad042SBarry Smith   PetscMPIInt    size,rank;
400273d9f13SBarry Smith 
401273d9f13SBarry Smith   PetscFunctionBegin;
402273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
403273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
404273d9f13SBarry Smith 
405b0a32e0cSBarry Smith   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
406b0a32e0cSBarry Smith   B->data             = (void*)b;
407273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
408273d9f13SBarry Smith   B->factor           = 0;
409273d9f13SBarry Smith   B->mapping          = 0;
410273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
411273d9f13SBarry Smith 
412ae77e121SBarry Smith   ierr = PetscSplitOwnership(B->comm,&B->m,&B->M);CHKERRQ(ierr);
413ae77e121SBarry Smith   B->n = B->N = PetscMax(B->N,B->n);CHKERRQ(ierr);
414273d9f13SBarry Smith 
415273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
416273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
4178a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
418273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
4198a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
420273d9f13SBarry Smith 
421b24ad042SBarry Smith   ierr = PetscMalloc((size+1)*sizeof(PetscInt),&b->rowners);CHKERRQ(ierr);
42252e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));CHKERRQ(ierr);
423273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
424273d9f13SBarry Smith   b->rowners[0] = 0;
425273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
426273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
427273d9f13SBarry Smith   }
428273d9f13SBarry Smith   b->rstart = b->rowners[rank];
429273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
430a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C",
431a23d5eceSKris Buschelman                                     "MatMPIAdjSetPreallocation_MPIAdj",
432a23d5eceSKris Buschelman                                      MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
433273d9f13SBarry Smith   PetscFunctionReturn(0);
434273d9f13SBarry Smith }
435273d9f13SBarry Smith EXTERN_C_END
436273d9f13SBarry Smith 
4374a2ae208SSatish Balay #undef __FUNCT__
4384a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
439a23d5eceSKris Buschelman /*@C
440a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
441a23d5eceSKris Buschelman 
442a23d5eceSKris Buschelman    Collective on MPI_Comm
443a23d5eceSKris Buschelman 
444a23d5eceSKris Buschelman    Input Parameters:
445a23d5eceSKris Buschelman +  A - the matrix
446a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
447a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
448a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
449a23d5eceSKris Buschelman -  values - [optional] edge weights
450a23d5eceSKris Buschelman 
451a23d5eceSKris Buschelman    Level: intermediate
452a23d5eceSKris Buschelman 
453a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
454a23d5eceSKris Buschelman @*/
455*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
456273d9f13SBarry Smith {
457b24ad042SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*,PetscInt*);
458273d9f13SBarry Smith 
459273d9f13SBarry Smith   PetscFunctionBegin;
460a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
461a23d5eceSKris Buschelman   if (f) {
462a23d5eceSKris Buschelman     ierr = (*f)(B,i,j,values);CHKERRQ(ierr);
463273d9f13SBarry Smith   }
464273d9f13SBarry Smith   PetscFunctionReturn(0);
465273d9f13SBarry Smith }
466273d9f13SBarry Smith 
4674a2ae208SSatish Balay #undef __FUNCT__
4684a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
469b97cf49bSBarry Smith /*@C
4703eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
471b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
472b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
473b97cf49bSBarry Smith 
474ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
475ef5ee4d1SLois Curfman McInnes 
476b97cf49bSBarry Smith    Input Parameters:
477c2e958feSBarry Smith +  comm - MPI communicator
4780752156aSBarry Smith .  m - number of local rows
479b97cf49bSBarry Smith .  n - number of columns
480b97cf49bSBarry Smith .  i - the indices into j for the start of each row
481329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
482ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
483329f5518SBarry Smith -  values -[optional] edge weights
484b97cf49bSBarry Smith 
485b97cf49bSBarry Smith    Output Parameter:
486b97cf49bSBarry Smith .  A - the matrix
487b97cf49bSBarry Smith 
48815091d37SBarry Smith    Level: intermediate
48915091d37SBarry Smith 
4904bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
4914bc6d8c0SBarry Smith    MatSetValues().
492c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
493ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
4941198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
495327e1a73SBarry Smith    Should not include the matrix diagonals.
496b97cf49bSBarry Smith 
497e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
498ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
499ededeb1aSvictorle 
500ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
501b97cf49bSBarry Smith 
502e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
503b97cf49bSBarry Smith @*/
504*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
505b97cf49bSBarry Smith {
506dfbe8321SBarry Smith   PetscErrorCode ierr;
507b97cf49bSBarry Smith 
508433994e6SBarry Smith   PetscFunctionBegin;
509273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
510273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
511273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
5123a40ed3dSBarry Smith   PetscFunctionReturn(0);
513b97cf49bSBarry Smith }
514b97cf49bSBarry Smith 
515273d9f13SBarry Smith EXTERN_C_BEGIN
5164a2ae208SSatish Balay #undef __FUNCT__
5174a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj"
518*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvertTo_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
519b3cc6726SBarry Smith {
520676c34cdSKris Buschelman   Mat               B;
5216849ba73SBarry Smith   PetscErrorCode    ierr;
522b24ad042SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
523b24ad042SBarry Smith   const PetscInt    *rj;
524b3cc6726SBarry Smith   const PetscScalar *ra;
525c2e958feSBarry Smith   MPI_Comm          comm;
526c2e958feSBarry Smith 
527c2e958feSBarry Smith   PetscFunctionBegin;
528c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
529c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
530c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
531c2e958feSBarry Smith 
532c2e958feSBarry Smith   /* count the number of nonzeros per row */
533c2e958feSBarry Smith   for (i=0; i<m; i++) {
534c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
535c2e958feSBarry Smith     for (j=0; j<len; j++) {
536c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
537c2e958feSBarry Smith     }
538c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
539c2e958feSBarry Smith     nzeros += len;
540c2e958feSBarry Smith   }
541c2e958feSBarry Smith 
542c2e958feSBarry Smith   /* malloc space for nonzeros */
543b24ad042SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr);
544b24ad042SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
545b24ad042SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr);
546c2e958feSBarry Smith 
547c2e958feSBarry Smith   nzeros = 0;
548c2e958feSBarry Smith   ia[0]  = 0;
549c2e958feSBarry Smith   for (i=0; i<m; i++) {
550c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
551c2e958feSBarry Smith     cnt     = 0;
552c2e958feSBarry Smith     for (j=0; j<len; j++) {
553c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
554b24ad042SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
555c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
556c2e958feSBarry Smith       }
557c2e958feSBarry Smith     }
558c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
559c2e958feSBarry Smith     nzeros += cnt;
560c2e958feSBarry Smith     ia[i+1] = nzeros;
561c2e958feSBarry Smith   }
562c2e958feSBarry Smith 
563c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
564f204ca49SKris Buschelman   ierr = MatCreate(comm,m,N,PETSC_DETERMINE,N,&B);CHKERRQ(ierr);
565f204ca49SKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
566f204ca49SKris Buschelman   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
567c2e958feSBarry Smith 
568ceb03754SKris Buschelman   if (reuse == MAT_REUSE_MATRIX) {
5698ab5b326SKris Buschelman     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
570c3d102feSKris Buschelman   } else {
571676c34cdSKris Buschelman     *newmat = B;
572c3d102feSKris Buschelman   }
573c2e958feSBarry Smith   PetscFunctionReturn(0);
574c2e958feSBarry Smith }
575273d9f13SBarry Smith EXTERN_C_END
576433994e6SBarry Smith 
577433994e6SBarry Smith 
578433994e6SBarry Smith 
579433994e6SBarry Smith 
580433994e6SBarry Smith 
581