xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 32077d6df18512e94702c86e7d562481ed0cebd0)
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"
9b0a32e0cSBarry Smith int MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)
10b97cf49bSBarry Smith {
113eda8832SBarry Smith   Mat_MPIAdj        *a = (Mat_MPIAdj*)A->data;
12fb9695e5SSatish Balay   int               ierr,i,j,m = A->m;
13fb9695e5SSatish Balay   char              *name;
14ce1411ecSBarry Smith   PetscViewerFormat format;
15b97cf49bSBarry Smith 
16433994e6SBarry Smith   PetscFunctionBegin;
173a7fca6bSBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
18b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
19fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
203a40ed3dSBarry Smith     PetscFunctionReturn(0);
21fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
22ffac6cdbSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
230752156aSBarry Smith   } else {
24b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
25b97cf49bSBarry Smith     for (i=0; i<m; i++) {
26b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr);
27b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
28b0a32e0cSBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr);
290752156aSBarry Smith       }
30b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
31b97cf49bSBarry Smith     }
32b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
33b97cf49bSBarry Smith   }
34b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
353a40ed3dSBarry Smith   PetscFunctionReturn(0);
36b97cf49bSBarry Smith }
37b97cf49bSBarry Smith 
384a2ae208SSatish Balay #undef __FUNCT__
394a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj"
40b0a32e0cSBarry Smith int MatView_MPIAdj(Mat A,PetscViewer viewer)
41b97cf49bSBarry Smith {
42b97cf49bSBarry Smith   int        ierr;
43*32077d6dSBarry Smith   PetscTruth iascii;
44b97cf49bSBarry Smith 
45433994e6SBarry Smith   PetscFunctionBegin;
46*32077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
47*32077d6dSBarry Smith   if (iascii) {
483eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
495cd90555SBarry Smith   } else {
5029bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name);
51b97cf49bSBarry Smith   }
523a40ed3dSBarry Smith   PetscFunctionReturn(0);
53b97cf49bSBarry Smith }
54b97cf49bSBarry Smith 
554a2ae208SSatish Balay #undef __FUNCT__
564a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj"
573eda8832SBarry Smith int MatDestroy_MPIAdj(Mat mat)
58b97cf49bSBarry Smith {
593eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data;
6094d884c6SBarry Smith   int        ierr;
6194d884c6SBarry Smith 
6294d884c6SBarry Smith   PetscFunctionBegin;
63aa482453SBarry Smith #if defined(PETSC_USE_LOG)
64b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz);
65b97cf49bSBarry Smith #endif
66606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
670400557aSBarry Smith   if (a->freeaij) {
68606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
69606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
701198db28SBarry Smith     if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);}
710400557aSBarry Smith   }
720400557aSBarry Smith   ierr = PetscFree(a->rowners);CHKERRQ(ierr);
731198db28SBarry Smith   ierr = PetscFree(a);CHKERRQ(ierr);
743a40ed3dSBarry Smith   PetscFunctionReturn(0);
75b97cf49bSBarry Smith }
76b97cf49bSBarry Smith 
774a2ae208SSatish Balay #undef __FUNCT__
784a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
793eda8832SBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op)
80b97cf49bSBarry Smith {
813eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
82b97cf49bSBarry Smith 
83433994e6SBarry Smith   PetscFunctionBegin;
8412c028f9SKris Buschelman   switch (op) {
8577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
8612c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
879a4540c5SBarry Smith   case MAT_HERMITIAN:
88b97cf49bSBarry Smith     a->symmetric = PETSC_TRUE;
8912c028f9SKris Buschelman     break;
909a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
919a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
929a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
939a4540c5SBarry Smith     a->symmetric = PETSC_FALSE;
949a4540c5SBarry Smith     break;
959a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
969a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
979a4540c5SBarry Smith     break;
9812c028f9SKris Buschelman   default:
99b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n");
10012c028f9SKris Buschelman     break;
101b97cf49bSBarry Smith   }
1023a40ed3dSBarry Smith   PetscFunctionReturn(0);
103b97cf49bSBarry Smith }
104b97cf49bSBarry Smith 
105b97cf49bSBarry Smith 
106b97cf49bSBarry Smith /*
107b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
108b97cf49bSBarry Smith */
109b97cf49bSBarry Smith 
1104a2ae208SSatish Balay #undef __FUNCT__
1114a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
1123eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A)
113b97cf49bSBarry Smith {
1143eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
11582502324SSatish Balay   int        i,j,*diag,m = A->m,ierr;
116b97cf49bSBarry Smith 
117433994e6SBarry Smith   PetscFunctionBegin;
118b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
119b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
120273d9f13SBarry Smith   for (i=0; i<A->m; i++) {
121b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
122b97cf49bSBarry Smith       if (a->j[j] == i) {
123b97cf49bSBarry Smith         diag[i] = j;
124b97cf49bSBarry Smith         break;
125b97cf49bSBarry Smith       }
126b97cf49bSBarry Smith     }
127b97cf49bSBarry Smith   }
128b97cf49bSBarry Smith   a->diag = diag;
1293a40ed3dSBarry Smith   PetscFunctionReturn(0);
130b97cf49bSBarry Smith }
131b97cf49bSBarry Smith 
1324a2ae208SSatish Balay #undef __FUNCT__
1334a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj"
13487828ca2SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v)
135b97cf49bSBarry Smith {
1363eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj*)A->data;
137b97cf49bSBarry Smith   int        *itmp;
138b97cf49bSBarry Smith 
139433994e6SBarry Smith   PetscFunctionBegin;
1400752156aSBarry Smith   row -= a->rstart;
1410752156aSBarry Smith 
142273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
143b97cf49bSBarry Smith 
144b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
145b97cf49bSBarry Smith   if (v) *v = PETSC_NULL;
146b97cf49bSBarry Smith   if (idx) {
147b97cf49bSBarry Smith     itmp = a->j + a->i[row];
148b97cf49bSBarry Smith     if (*nz) {
149b97cf49bSBarry Smith       *idx = itmp;
150b97cf49bSBarry Smith     }
151b97cf49bSBarry Smith     else *idx = 0;
152b97cf49bSBarry Smith   }
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
154b97cf49bSBarry Smith }
155b97cf49bSBarry Smith 
1564a2ae208SSatish Balay #undef __FUNCT__
1574a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj"
15887828ca2SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v)
159b97cf49bSBarry Smith {
160433994e6SBarry Smith   PetscFunctionBegin;
1613a40ed3dSBarry Smith   PetscFunctionReturn(0);
162b97cf49bSBarry Smith }
163b97cf49bSBarry Smith 
1644a2ae208SSatish Balay #undef __FUNCT__
1654a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj"
1663eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs)
167b97cf49bSBarry Smith {
168433994e6SBarry Smith   PetscFunctionBegin;
169b97cf49bSBarry Smith   *bs = 1;
1703a40ed3dSBarry Smith   PetscFunctionReturn(0);
171b97cf49bSBarry Smith }
172b97cf49bSBarry Smith 
173b97cf49bSBarry Smith 
1744a2ae208SSatish Balay #undef __FUNCT__
1754a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj"
1763eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg)
177b97cf49bSBarry Smith {
1783eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data;
1790f5bd95cSBarry Smith   int         ierr;
1800f5bd95cSBarry Smith   PetscTruth  flag;
181b97cf49bSBarry Smith 
182433994e6SBarry Smith   PetscFunctionBegin;
183b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
184273d9f13SBarry Smith   if ((A->m != B->m) ||(a->nz != b->nz)) {
1850f5bd95cSBarry Smith     flag = PETSC_FALSE;
186b97cf49bSBarry Smith   }
187b97cf49bSBarry Smith 
188b97cf49bSBarry Smith   /* if the a->i are the same */
189273d9f13SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr);
190b97cf49bSBarry Smith 
191b97cf49bSBarry Smith   /* if a->j are the same */
1920f5bd95cSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr);
193b97cf49bSBarry Smith 
194ca161407SBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr);
1953a40ed3dSBarry Smith   PetscFunctionReturn(0);
196b97cf49bSBarry Smith }
197b97cf49bSBarry Smith 
1984a2ae208SSatish Balay #undef __FUNCT__
1994a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
2004e7234bfSBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int *ia[],int *ja[],PetscTruth *done)
2016cd91f64SBarry Smith {
202d42735a1SBarry Smith   int        ierr,size,i;
2036cd91f64SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
2046cd91f64SBarry Smith 
2056cd91f64SBarry Smith   PetscFunctionBegin;
2066cd91f64SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
2076cd91f64SBarry Smith   if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);}
2086cd91f64SBarry Smith   *m    = A->m;
2096cd91f64SBarry Smith   *ia   = a->i;
2106cd91f64SBarry Smith   *ja   = a->j;
2116cd91f64SBarry Smith   *done = PETSC_TRUE;
212d42735a1SBarry Smith   if (oshift) {
213d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
214d42735a1SBarry Smith       (*ja)[i]++;
215d42735a1SBarry Smith     }
216d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
217d42735a1SBarry Smith   }
218d42735a1SBarry Smith   PetscFunctionReturn(0);
219d42735a1SBarry Smith }
220d42735a1SBarry Smith 
2214a2ae208SSatish Balay #undef __FUNCT__
2224a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
2234e7234bfSBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int *ia[],int *ja[],PetscTruth *done)
224d42735a1SBarry Smith {
225d42735a1SBarry Smith   int        i;
226d42735a1SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
227d42735a1SBarry Smith 
228d42735a1SBarry Smith   PetscFunctionBegin;
229c5b3d447SBarry Smith   if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()");
230c5b3d447SBarry Smith   if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()");
231d42735a1SBarry Smith   if (oshift) {
232d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
233d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
234d42735a1SBarry Smith       (*ja)[i]--;
235d42735a1SBarry Smith     }
236d42735a1SBarry Smith   }
2376cd91f64SBarry Smith   PetscFunctionReturn(0);
2386cd91f64SBarry Smith }
239b97cf49bSBarry Smith 
240b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2410b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
2423eda8832SBarry Smith        MatGetRow_MPIAdj,
2433eda8832SBarry Smith        MatRestoreRow_MPIAdj,
244b97cf49bSBarry Smith        0,
24597304618SKris Buschelman /* 4*/ 0,
246b97cf49bSBarry Smith        0,
247b97cf49bSBarry Smith        0,
248b97cf49bSBarry Smith        0,
2490b3b1487SBarry Smith        0,
2500b3b1487SBarry Smith        0,
25197304618SKris Buschelman /*10*/ 0,
2520b3b1487SBarry Smith        0,
2530b3b1487SBarry Smith        0,
2540b3b1487SBarry Smith        0,
2550b3b1487SBarry Smith        0,
25697304618SKris Buschelman /*15*/ 0,
2573eda8832SBarry Smith        MatEqual_MPIAdj,
2580b3b1487SBarry Smith        0,
2590b3b1487SBarry Smith        0,
2600b3b1487SBarry Smith        0,
26197304618SKris Buschelman /*20*/ 0,
2620b3b1487SBarry Smith        0,
2630b3b1487SBarry Smith        0,
2643eda8832SBarry Smith        MatSetOption_MPIAdj,
2650b3b1487SBarry Smith        0,
26697304618SKris Buschelman /*25*/ 0,
2670b3b1487SBarry Smith        0,
2680b3b1487SBarry Smith        0,
2690b3b1487SBarry Smith        0,
2700b3b1487SBarry Smith        0,
27197304618SKris Buschelman /*30*/ 0,
2720b3b1487SBarry Smith        0,
273273d9f13SBarry Smith        0,
274273d9f13SBarry Smith        0,
2750b3b1487SBarry Smith        0,
27697304618SKris Buschelman /*35*/ 0,
2770b3b1487SBarry Smith        0,
2780b3b1487SBarry Smith        0,
2790b3b1487SBarry Smith        0,
2800b3b1487SBarry Smith        0,
28197304618SKris Buschelman /*40*/ 0,
2820b3b1487SBarry Smith        0,
2830b3b1487SBarry Smith        0,
2840b3b1487SBarry Smith        0,
2850b3b1487SBarry Smith        0,
28697304618SKris Buschelman /*45*/ 0,
2870b3b1487SBarry Smith        0,
2880b3b1487SBarry Smith        0,
2890b3b1487SBarry Smith        0,
2900b3b1487SBarry Smith        0,
29197304618SKris Buschelman /*50*/ MatGetBlockSize_MPIAdj,
2926cd91f64SBarry Smith        MatGetRowIJ_MPIAdj,
293d42735a1SBarry Smith        MatRestoreRowIJ_MPIAdj,
294b97cf49bSBarry Smith        0,
295b97cf49bSBarry Smith        0,
29697304618SKris Buschelman /*55*/ 0,
297b97cf49bSBarry Smith        0,
2980752156aSBarry Smith        0,
2990752156aSBarry Smith        0,
3000b3b1487SBarry Smith        0,
30197304618SKris Buschelman /*60*/ 0,
302b9b97703SBarry Smith        MatDestroy_MPIAdj,
303b9b97703SBarry Smith        MatView_MPIAdj,
30497304618SKris Buschelman        MatGetPetscMaps_Petsc,
30597304618SKris Buschelman        0,
30697304618SKris Buschelman /*65*/ 0,
30797304618SKris Buschelman        0,
30897304618SKris Buschelman        0,
30997304618SKris Buschelman        0,
31097304618SKris Buschelman        0,
31197304618SKris Buschelman /*70*/ 0,
31297304618SKris Buschelman        0,
31397304618SKris Buschelman        0,
31497304618SKris Buschelman        0,
31597304618SKris Buschelman        0,
31697304618SKris Buschelman /*75*/ 0,
31797304618SKris Buschelman        0,
31897304618SKris Buschelman        0,
31997304618SKris Buschelman        0,
32097304618SKris Buschelman        0,
32197304618SKris Buschelman /*80*/ 0,
32297304618SKris Buschelman        0,
32397304618SKris Buschelman        0,
32497304618SKris Buschelman        0,
32597304618SKris Buschelman /*85*/ 0
32697304618SKris Buschelman };
327b97cf49bSBarry Smith 
328a23d5eceSKris Buschelman EXTERN_C_BEGIN
329a23d5eceSKris Buschelman #undef __FUNCT__
330a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
331a23d5eceSKris Buschelman int MatMPIAdjSetPreallocation_MPIAdj(Mat B,int *i,int *j,int *values)
332a23d5eceSKris Buschelman {
333a23d5eceSKris Buschelman   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
334a23d5eceSKris Buschelman   int        ierr;
335a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g)
336a23d5eceSKris Buschelman   int        ii;
337a23d5eceSKris Buschelman #endif
338a23d5eceSKris Buschelman 
339a23d5eceSKris Buschelman   PetscFunctionBegin;
340a23d5eceSKris Buschelman   B->preallocated = PETSC_TRUE;
341a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g)
342a23d5eceSKris Buschelman   if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]);
343a23d5eceSKris Buschelman   for (ii=1; ii<B->m; ii++) {
344a23d5eceSKris Buschelman     if (i[ii] < 0 || i[ii] < i[ii-1]) {
345a23d5eceSKris Buschelman       SETERRQ4(1,"i[%d]=%d index is out of range: i[%d]=%d",ii,i[ii],ii-1,i[ii-1]);
346a23d5eceSKris Buschelman     }
347a23d5eceSKris Buschelman   }
348a23d5eceSKris Buschelman   for (ii=0; ii<i[B->m]; ii++) {
349a23d5eceSKris Buschelman     if (j[ii] < 0 || j[ii] >= B->N) {
350a23d5eceSKris Buschelman       SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]);
351a23d5eceSKris Buschelman     }
352a23d5eceSKris Buschelman   }
353a23d5eceSKris Buschelman #endif
354a23d5eceSKris Buschelman 
355a23d5eceSKris Buschelman   b->j      = j;
356a23d5eceSKris Buschelman   b->i      = i;
357a23d5eceSKris Buschelman   b->values = values;
358a23d5eceSKris Buschelman 
359a23d5eceSKris Buschelman   b->nz               = i[B->m];
360a23d5eceSKris Buschelman   b->diag             = 0;
361a23d5eceSKris Buschelman   b->symmetric        = PETSC_FALSE;
362a23d5eceSKris Buschelman   b->freeaij          = PETSC_TRUE;
363a23d5eceSKris Buschelman 
364a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
365a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
366a23d5eceSKris Buschelman   PetscFunctionReturn(0);
367a23d5eceSKris Buschelman }
368a23d5eceSKris Buschelman EXTERN_C_END
369b97cf49bSBarry Smith 
3700bad9183SKris Buschelman /*MC
371fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
3720bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
3730bad9183SKris Buschelman 
3740bad9183SKris Buschelman   Level: beginner
3750bad9183SKris Buschelman 
3760bad9183SKris Buschelman .seealso: MatCreateMPIAdj
3770bad9183SKris Buschelman M*/
3780bad9183SKris Buschelman 
379273d9f13SBarry Smith EXTERN_C_BEGIN
3804a2ae208SSatish Balay #undef __FUNCT__
3814a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
382273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B)
383273d9f13SBarry Smith {
384273d9f13SBarry Smith   Mat_MPIAdj *b;
385273d9f13SBarry Smith   int        ii,ierr,size,rank;
386273d9f13SBarry Smith 
387273d9f13SBarry Smith   PetscFunctionBegin;
388273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
389273d9f13SBarry Smith   ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr);
390273d9f13SBarry Smith 
391b0a32e0cSBarry Smith   ierr                = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr);
392b0a32e0cSBarry Smith   B->data             = (void*)b;
393273d9f13SBarry Smith   ierr                = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr);
394273d9f13SBarry Smith   ierr                = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
395273d9f13SBarry Smith   B->factor           = 0;
396273d9f13SBarry Smith   B->lupivotthreshold = 1.0;
397273d9f13SBarry Smith   B->mapping          = 0;
398273d9f13SBarry Smith   B->assembled        = PETSC_FALSE;
399273d9f13SBarry Smith 
400273d9f13SBarry Smith   ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr);
401273d9f13SBarry Smith   B->n = B->N;
402273d9f13SBarry Smith 
403273d9f13SBarry Smith   /* the information in the maps duplicates the information computed below, eventually
404273d9f13SBarry Smith      we should remove the duplicate information that is not contained in the maps */
4058a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
406273d9f13SBarry Smith   /* we don't know the "local columns" so just use the row information :-(*/
4078a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr);
408273d9f13SBarry Smith 
409b0a32e0cSBarry Smith   ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr);
410b0a32e0cSBarry Smith   PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));
411273d9f13SBarry Smith   ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr);
412273d9f13SBarry Smith   b->rowners[0] = 0;
413273d9f13SBarry Smith   for (ii=2; ii<=size; ii++) {
414273d9f13SBarry Smith     b->rowners[ii] += b->rowners[ii-1];
415273d9f13SBarry Smith   }
416273d9f13SBarry Smith   b->rstart = b->rowners[rank];
417273d9f13SBarry Smith   b->rend   = b->rowners[rank+1];
418a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C",
419a23d5eceSKris Buschelman                                     "MatMPIAdjSetPreallocation_MPIAdj",
420a23d5eceSKris Buschelman                                      MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
421273d9f13SBarry Smith   PetscFunctionReturn(0);
422273d9f13SBarry Smith }
423273d9f13SBarry Smith EXTERN_C_END
424273d9f13SBarry Smith 
4254a2ae208SSatish Balay #undef __FUNCT__
4264a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
427a23d5eceSKris Buschelman /*@C
428a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
429a23d5eceSKris Buschelman 
430a23d5eceSKris Buschelman    Collective on MPI_Comm
431a23d5eceSKris Buschelman 
432a23d5eceSKris Buschelman    Input Parameters:
433a23d5eceSKris Buschelman +  A - the matrix
434a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
435a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
436a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
437a23d5eceSKris Buschelman -  values - [optional] edge weights
438a23d5eceSKris Buschelman 
439a23d5eceSKris Buschelman    Level: intermediate
440a23d5eceSKris Buschelman 
441a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
442a23d5eceSKris Buschelman @*/
443273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values)
444273d9f13SBarry Smith {
445a23d5eceSKris Buschelman   int ierr,(*f)(Mat,int*,int*,int*);
446273d9f13SBarry Smith 
447273d9f13SBarry Smith   PetscFunctionBegin;
448a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
449a23d5eceSKris Buschelman   if (f) {
450a23d5eceSKris Buschelman     ierr = (*f)(B,i,j,values);CHKERRQ(ierr);
451273d9f13SBarry Smith   }
452273d9f13SBarry Smith   PetscFunctionReturn(0);
453273d9f13SBarry Smith }
454273d9f13SBarry Smith 
4554a2ae208SSatish Balay #undef __FUNCT__
4564a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
457b97cf49bSBarry Smith /*@C
4583eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
459b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
460b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
461b97cf49bSBarry Smith 
462ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
463ef5ee4d1SLois Curfman McInnes 
464b97cf49bSBarry Smith    Input Parameters:
465c2e958feSBarry Smith +  comm - MPI communicator
4660752156aSBarry Smith .  m - number of local rows
467b97cf49bSBarry Smith .  n - number of columns
468b97cf49bSBarry Smith .  i - the indices into j for the start of each row
469329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
470ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
471329f5518SBarry Smith -  values -[optional] edge weights
472b97cf49bSBarry Smith 
473b97cf49bSBarry Smith    Output Parameter:
474b97cf49bSBarry Smith .  A - the matrix
475b97cf49bSBarry Smith 
47615091d37SBarry Smith    Level: intermediate
47715091d37SBarry Smith 
4784bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
4794bc6d8c0SBarry Smith    MatSetValues().
480c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
481ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
4821198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
483327e1a73SBarry Smith    Should not include the matrix diagonals.
484b97cf49bSBarry Smith 
485e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
486ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
487ededeb1aSvictorle 
488ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
489b97cf49bSBarry Smith 
490e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
491b97cf49bSBarry Smith @*/
4923eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A)
493b97cf49bSBarry Smith {
494273d9f13SBarry Smith   int        ierr;
495b97cf49bSBarry Smith 
496433994e6SBarry Smith   PetscFunctionBegin;
497273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr);
498273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
499273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
5003a40ed3dSBarry Smith   PetscFunctionReturn(0);
501b97cf49bSBarry Smith }
502b97cf49bSBarry Smith 
503273d9f13SBarry Smith EXTERN_C_BEGIN
5044a2ae208SSatish Balay #undef __FUNCT__
5054a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj"
506676c34cdSKris Buschelman int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *newmat) {
507676c34cdSKris Buschelman   Mat          B;
5084c49b128SBarry Smith   int          i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a;
50987828ca2SBarry Smith   PetscScalar  *ra;
510c2e958feSBarry Smith   MPI_Comm     comm;
511c2e958feSBarry Smith 
512c2e958feSBarry Smith   PetscFunctionBegin;
513c2e958feSBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
514c2e958feSBarry Smith   ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr);
515c2e958feSBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr);
516c2e958feSBarry Smith 
517c2e958feSBarry Smith   /* count the number of nonzeros per row */
518c2e958feSBarry Smith   for (i=0; i<m; i++) {
519c2e958feSBarry Smith     ierr   = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
520c2e958feSBarry Smith     for (j=0; j<len; j++) {
521c2e958feSBarry Smith       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
522c2e958feSBarry Smith     }
523c2e958feSBarry Smith     ierr   = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr);
524c2e958feSBarry Smith     nzeros += len;
525c2e958feSBarry Smith   }
526c2e958feSBarry Smith 
527c2e958feSBarry Smith   /* malloc space for nonzeros */
528b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr);
529b0a32e0cSBarry Smith   ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr);
530b0a32e0cSBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr);
531c2e958feSBarry Smith 
532c2e958feSBarry Smith   nzeros = 0;
533c2e958feSBarry Smith   ia[0]  = 0;
534c2e958feSBarry Smith   for (i=0; i<m; i++) {
535c2e958feSBarry Smith     ierr    = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
536c2e958feSBarry Smith     cnt     = 0;
537c2e958feSBarry Smith     for (j=0; j<len; j++) {
538c2e958feSBarry Smith       if (rj[j] != i+rstart) { /* if not diagonal */
5394c49b128SBarry Smith         a[nzeros+cnt]    = (int) PetscAbsScalar(ra[j]);
540c2e958feSBarry Smith         ja[nzeros+cnt++] = rj[j];
541c2e958feSBarry Smith       }
542c2e958feSBarry Smith     }
543c2e958feSBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
544c2e958feSBarry Smith     nzeros += cnt;
545c2e958feSBarry Smith     ia[i+1] = nzeros;
546c2e958feSBarry Smith   }
547c2e958feSBarry Smith 
548c2e958feSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
549f204ca49SKris Buschelman   ierr = MatCreate(comm,m,N,PETSC_DETERMINE,N,&B);CHKERRQ(ierr);
550f204ca49SKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
551f204ca49SKris Buschelman   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
552c2e958feSBarry Smith 
553676c34cdSKris Buschelman   /* Fake support for "inplace" convert. */
554676c34cdSKris Buschelman   if (*newmat == A) {
555676c34cdSKris Buschelman     ierr = MatDestroy(A);CHKERRQ(ierr);
556676c34cdSKris Buschelman   }
557676c34cdSKris Buschelman   *newmat = B;
558c2e958feSBarry Smith   PetscFunctionReturn(0);
559c2e958feSBarry Smith }
560273d9f13SBarry Smith EXTERN_C_END
561433994e6SBarry Smith 
562433994e6SBarry Smith 
563433994e6SBarry Smith 
564433994e6SBarry Smith 
565433994e6SBarry Smith 
566