xref: /petsc/src/mat/impls/is/matis.c (revision b4319ba4cb3cfcf3d1795d5e654729cd063fe205)
1*b4319ba4SBarry Smith /*$Id: is.c,v 1.15 2001/08/06 21:15:46 bsmith Exp $*/
2*b4319ba4SBarry Smith /*
3*b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
4*b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
5*b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
6*b4319ba4SBarry Smith    product is handled "implicitly".
7*b4319ba4SBarry Smith 
8*b4319ba4SBarry Smith      We provide:
9*b4319ba4SBarry Smith          MatMult()
10*b4319ba4SBarry Smith 
11*b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
12*b4319ba4SBarry Smith 
13*b4319ba4SBarry Smith */
14*b4319ba4SBarry Smith 
15*b4319ba4SBarry Smith #include "src/mat/impls/is/matis.h"      /*I "petscmat.h" I*/
16*b4319ba4SBarry Smith 
17*b4319ba4SBarry Smith #undef __FUNCT__
18*b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
19*b4319ba4SBarry Smith int MatDestroy_IS(Mat A)
20*b4319ba4SBarry Smith {
21*b4319ba4SBarry Smith   int    ierr;
22*b4319ba4SBarry Smith   Mat_IS *b = (Mat_IS*)A->data;
23*b4319ba4SBarry Smith 
24*b4319ba4SBarry Smith   PetscFunctionBegin;
25*b4319ba4SBarry Smith   if (b->A) {
26*b4319ba4SBarry Smith     ierr = MatDestroy(b->A);CHKERRQ(ierr);
27*b4319ba4SBarry Smith   }
28*b4319ba4SBarry Smith   if (b->ctx) {
29*b4319ba4SBarry Smith     ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr);
30*b4319ba4SBarry Smith   }
31*b4319ba4SBarry Smith   if (b->x) {
32*b4319ba4SBarry Smith     ierr = VecDestroy(b->x);CHKERRQ(ierr);
33*b4319ba4SBarry Smith   }
34*b4319ba4SBarry Smith   if (b->y) {
35*b4319ba4SBarry Smith     ierr = VecDestroy(b->y);CHKERRQ(ierr);
36*b4319ba4SBarry Smith   }
37*b4319ba4SBarry Smith   if (b->mapping) {
38*b4319ba4SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr);
39*b4319ba4SBarry Smith   }
40*b4319ba4SBarry Smith   ierr = PetscFree(b);CHKERRQ(ierr);
41*b4319ba4SBarry Smith   PetscFunctionReturn(0);
42*b4319ba4SBarry Smith }
43*b4319ba4SBarry Smith 
44*b4319ba4SBarry Smith #undef __FUNCT__
45*b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
46*b4319ba4SBarry Smith int MatMult_IS(Mat A,Vec x,Vec y)
47*b4319ba4SBarry Smith {
48*b4319ba4SBarry Smith   int         ierr;
49*b4319ba4SBarry Smith   Mat_IS      *is = (Mat_IS*)A->data;
50*b4319ba4SBarry Smith   PetscScalar zero = 0.0;
51*b4319ba4SBarry Smith 
52*b4319ba4SBarry Smith   PetscFunctionBegin;
53*b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
54*b4319ba4SBarry Smith   ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
55*b4319ba4SBarry Smith   ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
56*b4319ba4SBarry Smith 
57*b4319ba4SBarry Smith   /* multiply the local matrix */
58*b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
59*b4319ba4SBarry Smith 
60*b4319ba4SBarry Smith   /* scatter product back into global memory */
61*b4319ba4SBarry Smith   ierr = VecSet(&zero,y);CHKERRQ(ierr);
62*b4319ba4SBarry Smith   ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
63*b4319ba4SBarry Smith   ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
64*b4319ba4SBarry Smith 
65*b4319ba4SBarry Smith   PetscFunctionReturn(0);
66*b4319ba4SBarry Smith }
67*b4319ba4SBarry Smith 
68*b4319ba4SBarry Smith #undef __FUNCT__
69*b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
70*b4319ba4SBarry Smith int MatView_IS(Mat A,PetscViewer viewer)
71*b4319ba4SBarry Smith {
72*b4319ba4SBarry Smith   Mat_IS      *a = (Mat_IS*)A->data;
73*b4319ba4SBarry Smith   int         ierr;
74*b4319ba4SBarry Smith   PetscViewer sviewer;
75*b4319ba4SBarry Smith 
76*b4319ba4SBarry Smith   PetscFunctionBegin;
77*b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
78*b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
79*b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
80*b4319ba4SBarry Smith   PetscFunctionReturn(0);
81*b4319ba4SBarry Smith }
82*b4319ba4SBarry Smith 
83*b4319ba4SBarry Smith #undef __FUNCT__
84*b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
85*b4319ba4SBarry Smith int MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping)
86*b4319ba4SBarry Smith {
87*b4319ba4SBarry Smith   int    ierr,n;
88*b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)A->data;
89*b4319ba4SBarry Smith   IS     from,to;
90*b4319ba4SBarry Smith   Vec    global;
91*b4319ba4SBarry Smith 
92*b4319ba4SBarry Smith   PetscFunctionBegin;
93*b4319ba4SBarry Smith   is->mapping = mapping;
94*b4319ba4SBarry Smith   ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr);
95*b4319ba4SBarry Smith 
96*b4319ba4SBarry Smith   /* Create the local matrix A */
97*b4319ba4SBarry Smith   ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr);
98*b4319ba4SBarry Smith   ierr = MatCreate(PETSC_COMM_SELF,n,n,n,n,&is->A);CHKERRQ(ierr);
99*b4319ba4SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");CHKERRQ(ierr);
100*b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
101*b4319ba4SBarry Smith 
102*b4319ba4SBarry Smith   /* Create the local work vectors */
103*b4319ba4SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr);
104*b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
105*b4319ba4SBarry Smith 
106*b4319ba4SBarry Smith   /* setup the global to local scatter */
107*b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
108*b4319ba4SBarry Smith   ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr);
109*b4319ba4SBarry Smith   ierr = VecCreateMPI(A->comm,A->n,A->N,&global);CHKERRQ(ierr);
110*b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
111*b4319ba4SBarry Smith   ierr = VecDestroy(global);CHKERRQ(ierr);
112*b4319ba4SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
113*b4319ba4SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
114*b4319ba4SBarry Smith   PetscFunctionReturn(0);
115*b4319ba4SBarry Smith }
116*b4319ba4SBarry Smith 
117*b4319ba4SBarry Smith 
118*b4319ba4SBarry Smith #undef __FUNCT__
119*b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
120*b4319ba4SBarry Smith int MatSetValuesLocal_IS(Mat A,int m,const int *rows, int n,const int *cols,const PetscScalar *values,InsertMode addv)
121*b4319ba4SBarry Smith {
122*b4319ba4SBarry Smith   int    ierr;
123*b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)A->data;
124*b4319ba4SBarry Smith 
125*b4319ba4SBarry Smith   PetscFunctionBegin;
126*b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
127*b4319ba4SBarry Smith   PetscFunctionReturn(0);
128*b4319ba4SBarry Smith }
129*b4319ba4SBarry Smith 
130*b4319ba4SBarry Smith #undef __FUNCT__
131*b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
132*b4319ba4SBarry Smith int MatZeroRowsLocal_IS(Mat A,IS isrows,const PetscScalar *diag)
133*b4319ba4SBarry Smith {
134*b4319ba4SBarry Smith   Mat_IS      *is = (Mat_IS*)A->data;
135*b4319ba4SBarry Smith   int         ierr,i,n,*rows;
136*b4319ba4SBarry Smith   PetscScalar *array;
137*b4319ba4SBarry Smith 
138*b4319ba4SBarry Smith   PetscFunctionBegin;
139*b4319ba4SBarry Smith 
140*b4319ba4SBarry Smith   {
141*b4319ba4SBarry Smith     /*
142*b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
143*b4319ba4SBarry Smith        work properly in the interface nodes.
144*b4319ba4SBarry Smith     */
145*b4319ba4SBarry Smith     Vec         counter;
146*b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
147*b4319ba4SBarry Smith     ierr = VecCreateMPI(A->comm,A->n,A->N,&counter);CHKERRQ(ierr);
148*b4319ba4SBarry Smith     ierr = VecSet(&zero,counter);CHKERRQ(ierr);
149*b4319ba4SBarry Smith     ierr = VecSet(&one,is->x);CHKERRQ(ierr);
150*b4319ba4SBarry Smith     ierr = VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
151*b4319ba4SBarry Smith     ierr = VecScatterEnd  (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr);
152*b4319ba4SBarry Smith     ierr = VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
153*b4319ba4SBarry Smith     ierr = VecScatterEnd  (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr);
154*b4319ba4SBarry Smith     ierr = VecDestroy(counter);CHKERRQ(ierr);
155*b4319ba4SBarry Smith   }
156*b4319ba4SBarry Smith   ierr = ISGetLocalSize(isrows,&n);CHKERRQ(ierr);
157*b4319ba4SBarry Smith   if (n == 0) {
158*b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
159*b4319ba4SBarry Smith   } else {
160*b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
161*b4319ba4SBarry Smith     ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
162*b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
163*b4319ba4SBarry Smith     ierr = MatZeroRows(is->A,isrows,diag);CHKERRQ(ierr);
164*b4319ba4SBarry Smith     for (i=0; i<n; i++) {
165*b4319ba4SBarry Smith       ierr = MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
166*b4319ba4SBarry Smith     }
167*b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168*b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169*b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
170*b4319ba4SBarry Smith     ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
171*b4319ba4SBarry Smith   }
172*b4319ba4SBarry Smith 
173*b4319ba4SBarry Smith   PetscFunctionReturn(0);
174*b4319ba4SBarry Smith }
175*b4319ba4SBarry Smith 
176*b4319ba4SBarry Smith #undef __FUNCT__
177*b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
178*b4319ba4SBarry Smith int MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
179*b4319ba4SBarry Smith {
180*b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)A->data;
181*b4319ba4SBarry Smith   int    ierr;
182*b4319ba4SBarry Smith 
183*b4319ba4SBarry Smith   PetscFunctionBegin;
184*b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
185*b4319ba4SBarry Smith   PetscFunctionReturn(0);
186*b4319ba4SBarry Smith }
187*b4319ba4SBarry Smith 
188*b4319ba4SBarry Smith #undef __FUNCT__
189*b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
190*b4319ba4SBarry Smith int MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
191*b4319ba4SBarry Smith {
192*b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)A->data;
193*b4319ba4SBarry Smith   int    ierr;
194*b4319ba4SBarry Smith 
195*b4319ba4SBarry Smith   PetscFunctionBegin;
196*b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
197*b4319ba4SBarry Smith   PetscFunctionReturn(0);
198*b4319ba4SBarry Smith }
199*b4319ba4SBarry Smith 
200*b4319ba4SBarry Smith EXTERN_C_BEGIN
201*b4319ba4SBarry Smith #undef __FUNCT__
202*b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
203*b4319ba4SBarry Smith int MatISGetLocalMat_IS(Mat mat,Mat *local)
204*b4319ba4SBarry Smith {
205*b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS *)mat->data;
206*b4319ba4SBarry Smith 
207*b4319ba4SBarry Smith   PetscFunctionBegin;
208*b4319ba4SBarry Smith   *local = is->A;
209*b4319ba4SBarry Smith   PetscFunctionReturn(0);
210*b4319ba4SBarry Smith }
211*b4319ba4SBarry Smith EXTERN_C_END
212*b4319ba4SBarry Smith 
213*b4319ba4SBarry Smith #undef __FUNCT__
214*b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
215*b4319ba4SBarry Smith /*@
216*b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
217*b4319ba4SBarry Smith 
218*b4319ba4SBarry Smith   Input Parameter:
219*b4319ba4SBarry Smith .  mat - the matrix
220*b4319ba4SBarry Smith 
221*b4319ba4SBarry Smith   Output Parameter:
222*b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
223*b4319ba4SBarry Smith 
224*b4319ba4SBarry Smith   Level: advanced
225*b4319ba4SBarry Smith 
226*b4319ba4SBarry Smith   Notes:
227*b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
228*b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
229*b4319ba4SBarry Smith   of the MatSetValues() operation.
230*b4319ba4SBarry Smith 
231*b4319ba4SBarry Smith .seealso: MATIS
232*b4319ba4SBarry Smith @*/
233*b4319ba4SBarry Smith int MatISGetLocalMat(Mat mat,Mat *local)
234*b4319ba4SBarry Smith {
235*b4319ba4SBarry Smith   int ierr,(*f)(Mat,Mat *);
236*b4319ba4SBarry Smith 
237*b4319ba4SBarry Smith   PetscFunctionBegin;
238*b4319ba4SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
239*b4319ba4SBarry Smith   PetscValidPointer(local,2);
240*b4319ba4SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr);
241*b4319ba4SBarry Smith   if (f) {
242*b4319ba4SBarry Smith     ierr = (*f)(mat,local);CHKERRQ(ierr);
243*b4319ba4SBarry Smith   } else {
244*b4319ba4SBarry Smith     local = 0;
245*b4319ba4SBarry Smith   }
246*b4319ba4SBarry Smith   PetscFunctionReturn(0);
247*b4319ba4SBarry Smith }
248*b4319ba4SBarry Smith 
249*b4319ba4SBarry Smith /*MC
250*b4319ba4SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners.
251*b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
252*b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
253*b4319ba4SBarry Smith    product is handled "implicitly".
254*b4319ba4SBarry Smith 
255*b4319ba4SBarry Smith    Operations Provided:
256*b4319ba4SBarry Smith .  MatMult
257*b4319ba4SBarry Smith 
258*b4319ba4SBarry Smith    Options Database Keys:
259*b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
260*b4319ba4SBarry Smith 
261*b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
262*b4319ba4SBarry Smith 
263*b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
264*b4319ba4SBarry Smith 
265*b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
266*b4319ba4SBarry Smith           MatISGetLocalMat()
267*b4319ba4SBarry Smith 
268*b4319ba4SBarry Smith   Level: advanced
269*b4319ba4SBarry Smith 
270*b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping()
271*b4319ba4SBarry Smith 
272*b4319ba4SBarry Smith M*/
273*b4319ba4SBarry Smith 
274*b4319ba4SBarry Smith EXTERN_C_BEGIN
275*b4319ba4SBarry Smith #undef __FUNCT__
276*b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
277*b4319ba4SBarry Smith int MatCreate_IS(Mat A)
278*b4319ba4SBarry Smith {
279*b4319ba4SBarry Smith   int    ierr;
280*b4319ba4SBarry Smith   Mat_IS *b;
281*b4319ba4SBarry Smith 
282*b4319ba4SBarry Smith   PetscFunctionBegin;
283*b4319ba4SBarry Smith   ierr                = PetscNew(Mat_IS,&b);CHKERRQ(ierr);
284*b4319ba4SBarry Smith   A->data             = (void*)b;
285*b4319ba4SBarry Smith   ierr = PetscMemzero(b,sizeof(Mat_IS));CHKERRQ(ierr);
286*b4319ba4SBarry Smith   ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
287*b4319ba4SBarry Smith   A->factor           = 0;
288*b4319ba4SBarry Smith   A->mapping          = 0;
289*b4319ba4SBarry Smith 
290*b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
291*b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
292*b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
293*b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
294*b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
295*b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
296*b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
297*b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
298*b4319ba4SBarry Smith 
299*b4319ba4SBarry Smith   ierr = PetscSplitOwnership(A->comm,&A->m,&A->M);CHKERRQ(ierr);
300*b4319ba4SBarry Smith   ierr = PetscSplitOwnership(A->comm,&A->n,&A->N);CHKERRQ(ierr);
301*b4319ba4SBarry Smith   ierr = MPI_Scan(&A->m,&b->rend,1,MPI_INT,MPI_SUM,A->comm);CHKERRQ(ierr);
302*b4319ba4SBarry Smith   b->rstart = b->rend - A->m;
303*b4319ba4SBarry Smith 
304*b4319ba4SBarry Smith   b->A          = 0;
305*b4319ba4SBarry Smith   b->ctx        = 0;
306*b4319ba4SBarry Smith   b->x          = 0;
307*b4319ba4SBarry Smith   b->y          = 0;
308*b4319ba4SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr);
309*b4319ba4SBarry Smith 
310*b4319ba4SBarry Smith   PetscFunctionReturn(0);
311*b4319ba4SBarry Smith }
312*b4319ba4SBarry Smith EXTERN_C_END
313*b4319ba4SBarry Smith 
314*b4319ba4SBarry Smith 
315*b4319ba4SBarry Smith 
316*b4319ba4SBarry Smith 
317*b4319ba4SBarry Smith 
318*b4319ba4SBarry Smith 
319*b4319ba4SBarry Smith 
320*b4319ba4SBarry Smith 
321*b4319ba4SBarry Smith 
322*b4319ba4SBarry Smith 
323*b4319ba4SBarry Smith 
324*b4319ba4SBarry Smith 
325*b4319ba4SBarry Smith 
326