xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 2b1d87635c5979121742b260db309414b85f9c32)
1be1d678aSKris Buschelman 
2b97cf49bSBarry Smith /*
3b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
4b97cf49bSBarry Smith */
5e37d9405SJed Brown #include <../src/mat/impls/adj/mpi/mpiadj.h>    /*I "petscmat.h" I*/
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;
13d0f46423SBarry Smith   PetscInt          i,j,m = A->rmap->n;
142dcb1b2aSMatthew Knepley   const 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) {
23ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported");
240752156aSBarry Smith   } else {
25d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
267b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
27b97cf49bSBarry Smith     for (i=0; i<m; i++) {
28d0f46423SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr);
29b97cf49bSBarry Smith       for (j=a->i[i]; j<a->i[i+1]; j++) {
3077431f27SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr);
310752156aSBarry Smith       }
32b0a32e0cSBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
33b97cf49bSBarry Smith     }
34d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
35b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
367b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
377b23a99aSBarry Smith   }
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;
46ace3abfcSBarry Smith   PetscBool      iascii;
47b97cf49bSBarry Smith 
48433994e6SBarry Smith   PetscFunctionBegin;
49251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5032077d6dSBarry Smith   if (iascii) {
513eda8832SBarry Smith     ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr);
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)
65d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz);
66b97cf49bSBarry Smith #endif
6705b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
680400557aSBarry Smith   if (a->freeaij) {
6914196391SBarry Smith     if (a->freeaijwithfree) {
7014196391SBarry Smith       if (a->i) free(a->i);
7114196391SBarry Smith       if (a->j) free(a->j);
7214196391SBarry Smith     } else {
73606d414cSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
74606d414cSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
7505b42c5fSBarry Smith       ierr = PetscFree(a->values);CHKERRQ(ierr);
7614196391SBarry Smith     }
770400557aSBarry Smith   }
78*2b1d8763SJed Brown   ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
79bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
80dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
81bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr);
82bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr);
833a40ed3dSBarry Smith   PetscFunctionReturn(0);
84b97cf49bSBarry Smith }
85b97cf49bSBarry Smith 
864a2ae208SSatish Balay #undef __FUNCT__
874a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj"
88ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg)
89b97cf49bSBarry Smith {
903eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
9163ba0a88SBarry Smith   PetscErrorCode ierr;
92b97cf49bSBarry Smith 
93433994e6SBarry Smith   PetscFunctionBegin;
9412c028f9SKris Buschelman   switch (op) {
9577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
9612c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
979a4540c5SBarry Smith   case MAT_HERMITIAN:
984e0d8c25SBarry Smith     a->symmetric = flg;
999a4540c5SBarry Smith     break;
1009a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1019a4540c5SBarry Smith     break;
10212c028f9SKris Buschelman   default:
103290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
10412c028f9SKris Buschelman     break;
105b97cf49bSBarry Smith   }
1063a40ed3dSBarry Smith   PetscFunctionReturn(0);
107b97cf49bSBarry Smith }
108b97cf49bSBarry Smith 
109b97cf49bSBarry Smith 
110b97cf49bSBarry Smith /*
111b97cf49bSBarry Smith      Adds diagonal pointers to sparse matrix structure.
112b97cf49bSBarry Smith */
113b97cf49bSBarry Smith 
1144a2ae208SSatish Balay #undef __FUNCT__
1154a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj"
116dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A)
117b97cf49bSBarry Smith {
1183eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
1196849ba73SBarry Smith   PetscErrorCode ierr;
120d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n;
121b97cf49bSBarry Smith 
122433994e6SBarry Smith   PetscFunctionBegin;
12309f38230SBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
12409f38230SBarry Smith   ierr = PetscLogObjectMemory(A,m*sizeof(PetscInt));CHKERRQ(ierr);
125d0f46423SBarry Smith   for (i=0; i<A->rmap->n; i++) {
126b97cf49bSBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
127b97cf49bSBarry Smith       if (a->j[j] == i) {
12809f38230SBarry Smith         a->diag[i] = j;
129b97cf49bSBarry Smith         break;
130b97cf49bSBarry Smith       }
131b97cf49bSBarry Smith     }
132b97cf49bSBarry Smith   }
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;
141*2b1d8763SJed Brown   PetscErrorCode ierr;
142b97cf49bSBarry Smith 
143433994e6SBarry Smith   PetscFunctionBegin;
144d0f46423SBarry Smith   row -= A->rmap->rstart;
1450752156aSBarry Smith 
146e32f2f54SBarry Smith   if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
147b97cf49bSBarry Smith 
148b97cf49bSBarry Smith   *nz = a->i[row+1] - a->i[row];
149*2b1d8763SJed Brown   if (v) {
150*2b1d8763SJed Brown     PetscInt j;
151*2b1d8763SJed Brown     if (a->rowvalues_alloc < *nz) {
152*2b1d8763SJed Brown       ierr = PetscFree(a->rowvalues);CHKERRQ(ierr);
153*2b1d8763SJed Brown       a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz);
154*2b1d8763SJed Brown       ierr = PetscMalloc(a->rowvalues_alloc*sizeof(PetscScalar),&a->rowvalues);CHKERRQ(ierr);
155*2b1d8763SJed Brown     }
156*2b1d8763SJed Brown     for (j=0; j<*nz; j++) a->rowvalues[j] = a->values[a->i[row]+j];
157*2b1d8763SJed Brown     *v = (*nz) ? a->rowvalues : NULL;
158*2b1d8763SJed Brown   }
15992b6f2a0SJed Brown   if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
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"
173ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)
174b97cf49bSBarry Smith {
1753eda8832SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data;
176dfbe8321SBarry Smith   PetscErrorCode ierr;
177ace3abfcSBarry Smith   PetscBool      flag;
178b97cf49bSBarry Smith 
179433994e6SBarry Smith   PetscFunctionBegin;
180b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
181d0f46423SBarry Smith   if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) {
1820f5bd95cSBarry Smith     flag = PETSC_FALSE;
183b97cf49bSBarry Smith   }
184b97cf49bSBarry Smith 
185b97cf49bSBarry Smith   /* if the a->i are the same */
186d0f46423SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+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 
191ce94432eSBarry Smith   ierr = MPI_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1923a40ed3dSBarry Smith   PetscFunctionReturn(0);
193b97cf49bSBarry Smith }
194b97cf49bSBarry Smith 
1954a2ae208SSatish Balay #undef __FUNCT__
1964a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj"
1971a83f524SJed Brown PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
1986cd91f64SBarry Smith {
199b24ad042SBarry Smith   PetscInt   i;
2006cd91f64SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
2011a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
2026cd91f64SBarry Smith 
2036cd91f64SBarry Smith   PetscFunctionBegin;
204d0f46423SBarry Smith   *m    = A->rmap->n;
2056cd91f64SBarry Smith   *ia   = a->i;
2066cd91f64SBarry Smith   *ja   = a->j;
2076cd91f64SBarry Smith   *done = PETSC_TRUE;
208d42735a1SBarry Smith   if (oshift) {
209d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
210d42735a1SBarry Smith       (*ja)[i]++;
211d42735a1SBarry Smith     }
212d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]++;
213d42735a1SBarry Smith   }
214d42735a1SBarry Smith   PetscFunctionReturn(0);
215d42735a1SBarry Smith }
216d42735a1SBarry Smith 
2174a2ae208SSatish Balay #undef __FUNCT__
2184a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj"
2191a83f524SJed Brown PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
220d42735a1SBarry Smith {
221b24ad042SBarry Smith   PetscInt   i;
222d42735a1SBarry Smith   Mat_MPIAdj *a   = (Mat_MPIAdj*)A->data;
2231a83f524SJed Brown   PetscInt   **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
224d42735a1SBarry Smith 
225d42735a1SBarry Smith   PetscFunctionBegin;
226e32f2f54SBarry Smith   if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()");
227e32f2f54SBarry Smith   if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()");
228d42735a1SBarry Smith   if (oshift) {
229d42735a1SBarry Smith     for (i=0; i<=(*m); i++) (*ia)[i]--;
230d42735a1SBarry Smith     for (i=0; i<(*ia)[*m]; i++) {
231d42735a1SBarry Smith       (*ja)[i]--;
232d42735a1SBarry Smith     }
233d42735a1SBarry Smith   }
2346cd91f64SBarry Smith   PetscFunctionReturn(0);
2356cd91f64SBarry Smith }
236b97cf49bSBarry Smith 
23717667f90SBarry Smith #undef __FUNCT__
23817667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj"
23919fd82e9SBarry Smith PetscErrorCode  MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat)
24017667f90SBarry Smith {
24117667f90SBarry Smith   Mat               B;
24217667f90SBarry Smith   PetscErrorCode    ierr;
24317667f90SBarry Smith   PetscInt          i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a;
24417667f90SBarry Smith   const PetscInt    *rj;
24517667f90SBarry Smith   const PetscScalar *ra;
24617667f90SBarry Smith   MPI_Comm          comm;
24717667f90SBarry Smith 
24817667f90SBarry Smith   PetscFunctionBegin;
2490298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
2500298fd71SBarry Smith   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
2510298fd71SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
25217667f90SBarry Smith 
25317667f90SBarry Smith   /* count the number of nonzeros per row */
25417667f90SBarry Smith   for (i=0; i<m; i++) {
2550298fd71SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
2565ee9ba1cSJed Brown     for (j=0; j<len; j++) {
2575ee9ba1cSJed Brown       if (rj[j] == i+rstart) {len--; break;}    /* don't count diagonal */
2585ee9ba1cSJed Brown     }
25917667f90SBarry Smith     nzeros += len;
2600f2063bfSTobin Isaac     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr);
26117667f90SBarry Smith   }
26217667f90SBarry Smith 
26317667f90SBarry Smith   /* malloc space for nonzeros */
26417667f90SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr);
26517667f90SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr);
26617667f90SBarry Smith   ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr);
26717667f90SBarry Smith 
26817667f90SBarry Smith   nzeros = 0;
26917667f90SBarry Smith   ia[0]  = 0;
27017667f90SBarry Smith   for (i=0; i<m; i++) {
27117667f90SBarry Smith     ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
27217667f90SBarry Smith     cnt  = 0;
27317667f90SBarry Smith     for (j=0; j<len; j++) {
2745ee9ba1cSJed Brown       if (rj[j] != i+rstart) { /* if not diagonal */
27517667f90SBarry Smith         a[nzeros+cnt]    = (PetscInt) PetscAbsScalar(ra[j]);
27617667f90SBarry Smith         ja[nzeros+cnt++] = rj[j];
27717667f90SBarry Smith       }
2785ee9ba1cSJed Brown     }
27917667f90SBarry Smith     ierr    = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr);
28017667f90SBarry Smith     nzeros += cnt;
28117667f90SBarry Smith     ia[i+1] = nzeros;
28217667f90SBarry Smith   }
28317667f90SBarry Smith 
28417667f90SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
28517667f90SBarry Smith   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
28658ec18a5SLisandro Dalcin   ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
28717667f90SBarry Smith   ierr = MatSetType(B,type);CHKERRQ(ierr);
28817667f90SBarry Smith   ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr);
28917667f90SBarry Smith 
29017667f90SBarry Smith   if (reuse == MAT_REUSE_MATRIX) {
29117667f90SBarry Smith     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
29217667f90SBarry Smith   } else {
29317667f90SBarry Smith     *newmat = B;
29417667f90SBarry Smith   }
29517667f90SBarry Smith   PetscFunctionReturn(0);
29617667f90SBarry Smith }
29717667f90SBarry Smith 
298b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
2990b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0,
3003eda8832SBarry Smith                                        MatGetRow_MPIAdj,
3013eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
302b97cf49bSBarry Smith                                        0,
30397304618SKris Buschelman                                 /* 4*/ 0,
304b97cf49bSBarry Smith                                        0,
305b97cf49bSBarry Smith                                        0,
306b97cf49bSBarry Smith                                        0,
3070b3b1487SBarry Smith                                        0,
3080b3b1487SBarry Smith                                        0,
30997304618SKris Buschelman                                 /*10*/ 0,
3100b3b1487SBarry Smith                                        0,
3110b3b1487SBarry Smith                                        0,
3120b3b1487SBarry Smith                                        0,
3130b3b1487SBarry Smith                                        0,
31497304618SKris Buschelman                                 /*15*/ 0,
3153eda8832SBarry Smith                                        MatEqual_MPIAdj,
3160b3b1487SBarry Smith                                        0,
3170b3b1487SBarry Smith                                        0,
3180b3b1487SBarry Smith                                        0,
31997304618SKris Buschelman                                 /*20*/ 0,
3200b3b1487SBarry Smith                                        0,
3213eda8832SBarry Smith                                        MatSetOption_MPIAdj,
3220b3b1487SBarry Smith                                        0,
323d519adbfSMatthew Knepley                                 /*24*/ 0,
3240b3b1487SBarry Smith                                        0,
3250b3b1487SBarry Smith                                        0,
3260b3b1487SBarry Smith                                        0,
3270b3b1487SBarry Smith                                        0,
328d519adbfSMatthew Knepley                                 /*29*/ 0,
3290b3b1487SBarry Smith                                        0,
330273d9f13SBarry Smith                                        0,
331273d9f13SBarry Smith                                        0,
3320b3b1487SBarry Smith                                        0,
333d519adbfSMatthew Knepley                                 /*34*/ 0,
3340b3b1487SBarry Smith                                        0,
3350b3b1487SBarry Smith                                        0,
3360b3b1487SBarry Smith                                        0,
3370b3b1487SBarry Smith                                        0,
338d519adbfSMatthew Knepley                                 /*39*/ 0,
3390b3b1487SBarry Smith                                        0,
3400b3b1487SBarry Smith                                        0,
3410b3b1487SBarry Smith                                        0,
3420b3b1487SBarry Smith                                        0,
343d519adbfSMatthew Knepley                                 /*44*/ 0,
3440b3b1487SBarry Smith                                        0,
3450b3b1487SBarry Smith                                        0,
3460b3b1487SBarry Smith                                        0,
3470b3b1487SBarry Smith                                        0,
348d519adbfSMatthew Knepley                                 /*49*/ 0,
3496cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
350d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
351b97cf49bSBarry Smith                                        0,
352b97cf49bSBarry Smith                                        0,
353d519adbfSMatthew Knepley                                 /*54*/ 0,
354b97cf49bSBarry Smith                                        0,
3550752156aSBarry Smith                                        0,
3560752156aSBarry Smith                                        0,
3570b3b1487SBarry Smith                                        0,
358d519adbfSMatthew Knepley                                 /*59*/ 0,
359b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
360b9b97703SBarry Smith                                        MatView_MPIAdj,
36117667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
36297304618SKris Buschelman                                        0,
363d519adbfSMatthew Knepley                                 /*64*/ 0,
36497304618SKris Buschelman                                        0,
36597304618SKris Buschelman                                        0,
36697304618SKris Buschelman                                        0,
36797304618SKris Buschelman                                        0,
368d519adbfSMatthew Knepley                                 /*69*/ 0,
36997304618SKris Buschelman                                        0,
37097304618SKris Buschelman                                        0,
37197304618SKris Buschelman                                        0,
37297304618SKris Buschelman                                        0,
373d519adbfSMatthew Knepley                                 /*74*/ 0,
37497304618SKris Buschelman                                        0,
37597304618SKris Buschelman                                        0,
37697304618SKris Buschelman                                        0,
37797304618SKris Buschelman                                        0,
378d519adbfSMatthew Knepley                                 /*79*/ 0,
37997304618SKris Buschelman                                        0,
38097304618SKris Buschelman                                        0,
38197304618SKris Buschelman                                        0,
382865e5f61SKris Buschelman                                        0,
383d519adbfSMatthew Knepley                                 /*84*/ 0,
384865e5f61SKris Buschelman                                        0,
385865e5f61SKris Buschelman                                        0,
386865e5f61SKris Buschelman                                        0,
387865e5f61SKris Buschelman                                        0,
388d519adbfSMatthew Knepley                                 /*89*/ 0,
389865e5f61SKris Buschelman                                        0,
390865e5f61SKris Buschelman                                        0,
391865e5f61SKris Buschelman                                        0,
392865e5f61SKris Buschelman                                        0,
393d519adbfSMatthew Knepley                                 /*94*/ 0,
394865e5f61SKris Buschelman                                        0,
395865e5f61SKris Buschelman                                        0,
3963964eb88SJed Brown                                        0,
3973964eb88SJed Brown                                        0,
3983964eb88SJed Brown                                 /*99*/ 0,
3993964eb88SJed Brown                                        0,
4003964eb88SJed Brown                                        0,
4013964eb88SJed Brown                                        0,
4023964eb88SJed Brown                                        0,
4033964eb88SJed Brown                                /*104*/ 0,
4043964eb88SJed Brown                                        0,
4053964eb88SJed Brown                                        0,
4063964eb88SJed Brown                                        0,
4073964eb88SJed Brown                                        0,
4083964eb88SJed Brown                                /*109*/ 0,
4093964eb88SJed Brown                                        0,
4103964eb88SJed Brown                                        0,
4113964eb88SJed Brown                                        0,
4123964eb88SJed Brown                                        0,
4133964eb88SJed Brown                                /*114*/ 0,
4143964eb88SJed Brown                                        0,
4153964eb88SJed Brown                                        0,
4163964eb88SJed Brown                                        0,
4173964eb88SJed Brown                                        0,
4183964eb88SJed Brown                                /*119*/ 0,
4193964eb88SJed Brown                                        0,
4203964eb88SJed Brown                                        0,
4213964eb88SJed Brown                                        0,
4223964eb88SJed Brown                                        0,
4233964eb88SJed Brown                                /*124*/ 0,
4243964eb88SJed Brown                                        0,
4253964eb88SJed Brown                                        0,
4263964eb88SJed Brown                                        0,
4273964eb88SJed Brown                                        0,
4283964eb88SJed Brown                                /*129*/ 0,
4293964eb88SJed Brown                                        0,
4303964eb88SJed Brown                                        0,
4313964eb88SJed Brown                                        0,
4323964eb88SJed Brown                                        0,
4333964eb88SJed Brown                                /*134*/ 0,
4343964eb88SJed Brown                                        0,
4353964eb88SJed Brown                                        0,
4363964eb88SJed Brown                                        0,
4373964eb88SJed Brown                                        0,
4383964eb88SJed Brown                                /*139*/ 0,
4393964eb88SJed Brown                                        0
4403964eb88SJed Brown };
441b97cf49bSBarry Smith 
442a23d5eceSKris Buschelman #undef __FUNCT__
443a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj"
444f7a08781SBarry Smith static PetscErrorCode  MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
445a23d5eceSKris Buschelman {
446a23d5eceSKris Buschelman   Mat_MPIAdj     *b = (Mat_MPIAdj*)B->data;
447dfbe8321SBarry Smith   PetscErrorCode ierr;
4482515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
449b24ad042SBarry Smith   PetscInt ii;
450a23d5eceSKris Buschelman #endif
451a23d5eceSKris Buschelman 
452a23d5eceSKris Buschelman   PetscFunctionBegin;
45326283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
45426283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
45558ec18a5SLisandro Dalcin 
4562515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
457e32f2f54SBarry Smith   if (i[0] != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]);
458d0f46423SBarry Smith   for (ii=1; ii<B->rmap->n; ii++) {
459e02043d6SBarry Smith     if (i[ii] < 0 || i[ii] < i[ii-1]) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i[%D]=%D index is out of range: i[%D]=%D",ii,i[ii],ii-1,i[ii-1]);
460a23d5eceSKris Buschelman   }
461d0f46423SBarry Smith   for (ii=0; ii<i[B->rmap->n]; ii++) {
462e02043d6SBarry Smith     if (j[ii] < 0 || j[ii] >= B->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]);
463a23d5eceSKris Buschelman   }
464a23d5eceSKris Buschelman #endif
46558ec18a5SLisandro Dalcin   B->preallocated = PETSC_TRUE;
466a23d5eceSKris Buschelman 
467a23d5eceSKris Buschelman   b->j      = j;
468a23d5eceSKris Buschelman   b->i      = i;
469a23d5eceSKris Buschelman   b->values = values;
470a23d5eceSKris Buschelman 
471d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
472a23d5eceSKris Buschelman   b->diag      = 0;
473a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
474a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
475a23d5eceSKris Buschelman 
476a23d5eceSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
477a23d5eceSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
478a23d5eceSKris Buschelman   PetscFunctionReturn(0);
479a23d5eceSKris Buschelman }
480b97cf49bSBarry Smith 
4819a3665c6SJed Brown #undef __FUNCT__
4829a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj"
483f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B)
4849a3665c6SJed Brown {
4859a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj*)A->data;
4869a3665c6SJed Brown   PetscErrorCode ierr;
4879a3665c6SJed Brown   const PetscInt *ranges;
4889a3665c6SJed Brown   MPI_Comm       acomm,bcomm;
4899a3665c6SJed Brown   MPI_Group      agroup,bgroup;
4909a3665c6SJed Brown   PetscMPIInt    i,rank,size,nranks,*ranks;
4919a3665c6SJed Brown 
4929a3665c6SJed Brown   PetscFunctionBegin;
4930298fd71SBarry Smith   *B    = NULL;
494ce94432eSBarry Smith   ierr  = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr);
4959a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&size);CHKERRQ(ierr);
4969a3665c6SJed Brown   ierr  = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr);
4979a3665c6SJed Brown   ierr  = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr);
4989a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
4999a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) nranks++;
5009a3665c6SJed Brown   }
5019a3665c6SJed Brown   if (nranks == size) {         /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
5029a3665c6SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
5039a3665c6SJed Brown     *B   = A;
5049a3665c6SJed Brown     PetscFunctionReturn(0);
5059a3665c6SJed Brown   }
5069a3665c6SJed Brown 
5079a3665c6SJed Brown   ierr = PetscMalloc(nranks*sizeof(PetscMPIInt),&ranks);CHKERRQ(ierr);
5089a3665c6SJed Brown   for (i=0,nranks=0; i<size; i++) {
5099a3665c6SJed Brown     if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i;
5109a3665c6SJed Brown   }
5119a3665c6SJed Brown   ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr);
5129a3665c6SJed Brown   ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr);
5139a3665c6SJed Brown   ierr = PetscFree(ranks);CHKERRQ(ierr);
5149a3665c6SJed Brown   ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr);
5159a3665c6SJed Brown   ierr = MPI_Group_free(&agroup);CHKERRQ(ierr);
5169a3665c6SJed Brown   ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr);
5179a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
5189a3665c6SJed Brown     PetscInt   m,N;
5199a3665c6SJed Brown     Mat_MPIAdj *b;
5200298fd71SBarry Smith     ierr       = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
5210298fd71SBarry Smith     ierr       = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
5229a3665c6SJed Brown     ierr       = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr);
5239a3665c6SJed Brown     b          = (Mat_MPIAdj*)(*B)->data;
5249a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
5259a3665c6SJed Brown     ierr       = MPI_Comm_free(&bcomm);CHKERRQ(ierr);
5269a3665c6SJed Brown   }
5279a3665c6SJed Brown   PetscFunctionReturn(0);
5289a3665c6SJed Brown }
5299a3665c6SJed Brown 
5309a3665c6SJed Brown #undef __FUNCT__
5319a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat"
5329a3665c6SJed Brown /*@
5339a3665c6SJed Brown    MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows
5349a3665c6SJed Brown 
5359a3665c6SJed Brown    Collective
5369a3665c6SJed Brown 
5379a3665c6SJed Brown    Input Arguments:
5389a3665c6SJed Brown .  A - original MPIAdj matrix
5399a3665c6SJed Brown 
5409a3665c6SJed Brown    Output Arguments:
5410298fd71SBarry Smith .  B - matrix on subcommunicator, NULL on ranks that owned zero rows of A
5429a3665c6SJed Brown 
5439a3665c6SJed Brown    Level: developer
5449a3665c6SJed Brown 
5459a3665c6SJed Brown    Note:
5469a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
5479a3665c6SJed Brown 
5489a3665c6SJed Brown    The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed.
5499a3665c6SJed Brown 
5509a3665c6SJed Brown .seealso: MatCreateMPIAdj()
5519a3665c6SJed Brown @*/
5529a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B)
5539a3665c6SJed Brown {
5549a3665c6SJed Brown   PetscErrorCode ierr;
5559a3665c6SJed Brown 
5569a3665c6SJed Brown   PetscFunctionBegin;
5579a3665c6SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
5589a3665c6SJed Brown   ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr);
5599a3665c6SJed Brown   PetscFunctionReturn(0);
5609a3665c6SJed Brown }
5619a3665c6SJed Brown 
5620bad9183SKris Buschelman /*MC
563fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
5640bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
5650bad9183SKris Buschelman 
5660bad9183SKris Buschelman   Level: beginner
5670bad9183SKris Buschelman 
5680bad9183SKris Buschelman .seealso: MatCreateMPIAdj
5690bad9183SKris Buschelman M*/
5700bad9183SKris Buschelman 
5714a2ae208SSatish Balay #undef __FUNCT__
5724a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj"
5738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
574273d9f13SBarry Smith {
575273d9f13SBarry Smith   Mat_MPIAdj     *b;
5766849ba73SBarry Smith   PetscErrorCode ierr;
577273d9f13SBarry Smith 
578273d9f13SBarry Smith   PetscFunctionBegin;
57938f2d2fdSLisandro Dalcin   ierr         = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr);
580b0a32e0cSBarry Smith   B->data      = (void*)b;
581273d9f13SBarry Smith   ierr         = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
582273d9f13SBarry Smith   B->assembled = PETSC_FALSE;
583273d9f13SBarry Smith 
584bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr);
585bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr);
58617667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr);
587273d9f13SBarry Smith   PetscFunctionReturn(0);
588273d9f13SBarry Smith }
589273d9f13SBarry Smith 
5904a2ae208SSatish Balay #undef __FUNCT__
5914a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation"
592a23d5eceSKris Buschelman /*@C
593a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
594a23d5eceSKris Buschelman 
5953f9fe445SBarry Smith    Logically Collective on MPI_Comm
596a23d5eceSKris Buschelman 
597a23d5eceSKris Buschelman    Input Parameters:
598a23d5eceSKris Buschelman +  A - the matrix
599a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
600a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
601a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
602a23d5eceSKris Buschelman -  values - [optional] edge weights
603a23d5eceSKris Buschelman 
604a23d5eceSKris Buschelman    Level: intermediate
605a23d5eceSKris Buschelman 
606a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues()
607a23d5eceSKris Buschelman @*/
6087087cfbeSBarry Smith PetscErrorCode  MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values)
609273d9f13SBarry Smith {
6104ac538c5SBarry Smith   PetscErrorCode ierr;
611273d9f13SBarry Smith 
612273d9f13SBarry Smith   PetscFunctionBegin;
6134ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr);
614273d9f13SBarry Smith   PetscFunctionReturn(0);
615273d9f13SBarry Smith }
616273d9f13SBarry Smith 
6174a2ae208SSatish Balay #undef __FUNCT__
6184a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj"
619b97cf49bSBarry Smith /*@C
6203eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
621b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
622b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
623b97cf49bSBarry Smith 
624ef5ee4d1SLois Curfman McInnes    Collective on MPI_Comm
625ef5ee4d1SLois Curfman McInnes 
626b97cf49bSBarry Smith    Input Parameters:
627c2e958feSBarry Smith +  comm - MPI communicator
6280752156aSBarry Smith .  m - number of local rows
62958ec18a5SLisandro Dalcin .  N - number of global columns
630b97cf49bSBarry Smith .  i - the indices into j for the start of each row
631329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
632ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
633329f5518SBarry Smith -  values -[optional] edge weights
634b97cf49bSBarry Smith 
635b97cf49bSBarry Smith    Output Parameter:
636b97cf49bSBarry Smith .  A - the matrix
637b97cf49bSBarry Smith 
63815091d37SBarry Smith    Level: intermediate
63915091d37SBarry Smith 
6404bc6d8c0SBarry Smith    Notes: This matrix object does not support most matrix operations, include
6414bc6d8c0SBarry Smith    MatSetValues().
642c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
643ededeb1aSvictorle    when the matrix is destroyed; you must allocate them with PetscMalloc(). If you
6441198db28SBarry Smith     call from Fortran you need not create the arrays with PetscMalloc().
645327e1a73SBarry Smith    Should not include the matrix diagonals.
646b97cf49bSBarry Smith 
647e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
648ededeb1aSvictorle    to MatConvert, specifying a type of MATMPIADJ.
649ededeb1aSvictorle 
650ef5ee4d1SLois Curfman McInnes    Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC
651b97cf49bSBarry Smith 
652e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering()
653b97cf49bSBarry Smith @*/
6547087cfbeSBarry Smith PetscErrorCode  MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A)
655b97cf49bSBarry Smith {
656dfbe8321SBarry Smith   PetscErrorCode ierr;
657b97cf49bSBarry Smith 
658433994e6SBarry Smith   PetscFunctionBegin;
659f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
66058ec18a5SLisandro Dalcin   ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr);
661273d9f13SBarry Smith   ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr);
662273d9f13SBarry Smith   ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr);
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
664b97cf49bSBarry Smith }
665b97cf49bSBarry Smith 
666c2e958feSBarry Smith 
667433994e6SBarry Smith 
668433994e6SBarry Smith 
669433994e6SBarry Smith 
670433994e6SBarry Smith 
671433994e6SBarry Smith 
672