xref: /petsc/src/mat/utils/getcolv.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I   "petscmat.h"  I*/
20925cdddSBarry Smith 
30925cdddSBarry Smith /*@
482bf6240SBarry Smith   MatGetColumnVector - Gets the values from a given column of a matrix.
50925cdddSBarry Smith 
672257631SBarry Smith   Not Collective
7fee21e36SBarry Smith 
898a79cdbSBarry Smith   Input Parameters:
95ed6d96aSBarry Smith + A   - the matrix
105ed6d96aSBarry Smith . yy  - the vector
11b7e53d48SBarry Smith - col - the column requested (in global numbering)
1298a79cdbSBarry Smith 
1315091d37SBarry Smith   Level: advanced
1415091d37SBarry Smith 
159d006df2SBarry Smith   Notes:
1627430b45SBarry Smith   If a `MatType` does not implement the operation, each processor for which this is called
1711a5261eSBarry Smith   gets the values for its rows using `MatGetRow()`.
189d006df2SBarry Smith 
193d81755aSBarry Smith   The vector must have the same parallel row layout as the matrix.
203d81755aSBarry Smith 
2182bf6240SBarry Smith   Contributed by: Denis Vanderstraeten
220925cdddSBarry Smith 
231cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetRow()`, `MatGetDiagonal()`, `MatMult()`
240925cdddSBarry Smith @*/
MatGetColumnVector(Mat A,Vec yy,PetscInt col)25d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnVector(Mat A, Vec yy, PetscInt col)
26d71ae5a4SJacob Faibussowitsch {
278aa348c1SBarry Smith   PetscScalar       *y;
28b3cc6726SBarry Smith   const PetscScalar *v;
2938baddfdSBarry Smith   PetscInt           i, j, nz, N, Rs, Re, rs, re;
3038baddfdSBarry Smith   const PetscInt    *idx;
310925cdddSBarry Smith 
320925cdddSBarry Smith   PetscFunctionBegin;
330700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
340700a824SBarry Smith   PetscValidHeaderSpecific(yy, VEC_CLASSID, 2);
359a971ab9SStefano Zampini   PetscValidLogicalCollectiveInt(A, col, 3);
3608401ef6SPierre Jolivet   PetscCheck(col >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Requested negative column: %" PetscInt_FMT, col);
379566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &N));
3808401ef6SPierre Jolivet   PetscCheck(col < N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Requested column %" PetscInt_FMT " larger than number columns in matrix %" PetscInt_FMT, col, N);
399566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &Rs, &Re));
409566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(yy, &rs, &re));
41aed4548fSBarry Smith   PetscCheck(Rs == rs && Re == re, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Matrix %" PetscInt_FMT " %" PetscInt_FMT " does not have same ownership range (size) as vector %" PetscInt_FMT " %" PetscInt_FMT, Rs, Re, rs, re);
4282bf6240SBarry Smith 
43dbbe0bcdSBarry Smith   if (A->ops->getcolumnvector) PetscUseTypeMethod(A, getcolumnvector, yy, col);
44dbbe0bcdSBarry Smith   else {
459566063dSJacob Faibussowitsch     PetscCall(VecSet(yy, 0.0));
469566063dSJacob Faibussowitsch     PetscCall(VecGetArray(yy, &y));
479a971ab9SStefano Zampini     /* TODO for general matrices */
4882bf6240SBarry Smith     for (i = Rs; i < Re; i++) {
499566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A, i, &nz, &idx, &v));
5082bf6240SBarry Smith       if (nz && idx[0] <= col) {
5182bf6240SBarry Smith         /*
5282bf6240SBarry Smith           Should use faster search here
5382bf6240SBarry Smith         */
5482bf6240SBarry Smith         for (j = 0; j < nz; j++) {
5582bf6240SBarry Smith           if (idx[j] >= col) {
5682bf6240SBarry Smith             if (idx[j] == col) y[i - rs] = v[j];
5782bf6240SBarry Smith             break;
580925cdddSBarry Smith           }
590925cdddSBarry Smith         }
600925cdddSBarry Smith       }
619566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A, i, &nz, &idx, &v));
620925cdddSBarry Smith     }
639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(yy, &y));
648d0534beSBarry Smith   }
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
660925cdddSBarry Smith }
67242f1d38SBarry Smith 
6886819fdcSBarry Smith /*@
690716a85fSBarry Smith   MatGetColumnNorms - Gets the norms of each column of a sparse or dense matrix.
70242f1d38SBarry Smith 
71d8d19677SJose E. Roman   Input Parameters:
72242f1d38SBarry Smith + A    - the matrix
7311a5261eSBarry Smith - type - `NORM_2`, `NORM_1` or `NORM_INFINITY`
74242f1d38SBarry Smith 
75242f1d38SBarry Smith   Output Parameter:
76242f1d38SBarry Smith . norms - an array as large as the TOTAL number of columns in the matrix
77242f1d38SBarry Smith 
78f6680f47SSatish Balay   Level: intermediate
79f6680f47SSatish Balay 
8011a5261eSBarry Smith   Note:
8195452b02SPatrick Sanan   Each process has ALL the column norms after the call. Because of the way this is computed each process gets all the values,
8286819fdcSBarry Smith   if each process wants only some of the values it should extract the ones it wants from the array.
8386819fdcSBarry Smith 
841cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `NormType`, `MatNorm()`
8586819fdcSBarry Smith @*/
MatGetColumnNorms(Mat A,NormType type,PetscReal norms[])86d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnNorms(Mat A, NormType type, PetscReal norms[])
87d71ae5a4SJacob Faibussowitsch {
88857cbf51SRichard Tran Mills   /* NOTE: MatGetColumnNorms() could simply be a macro that calls MatGetColumnReductions().
89857cbf51SRichard Tran Mills    * I've kept this as a function because it allows slightly more in the way of error checking,
90857cbf51SRichard Tran Mills    * erroring out if MatGetColumnNorms() is not called with a valid NormType. */
91242f1d38SBarry Smith 
92242f1d38SBarry Smith   PetscFunctionBegin;
93*966bd95aSPierre Jolivet   PetscCheck(type == NORM_2 || type == NORM_1 || type == NORM_FROBENIUS || type == NORM_INFINITY || type == NORM_1_AND_2, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Unknown NormType");
949566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, type, norms));
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96a873a8cdSSam Reynolds }
97857cbf51SRichard Tran Mills 
98857cbf51SRichard Tran Mills /*@
99857cbf51SRichard Tran Mills   MatGetColumnSumsRealPart - Gets the sums of the real part of each column of a sparse or dense matrix.
100857cbf51SRichard Tran Mills 
101857cbf51SRichard Tran Mills   Input Parameter:
102857cbf51SRichard Tran Mills . A - the matrix
103857cbf51SRichard Tran Mills 
104857cbf51SRichard Tran Mills   Output Parameter:
105857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix
106857cbf51SRichard Tran Mills 
107857cbf51SRichard Tran Mills   Level: intermediate
108857cbf51SRichard Tran Mills 
10911a5261eSBarry Smith   Note:
110857cbf51SRichard Tran Mills   Each process has ALL the column sums after the call. Because of the way this is computed each process gets all the values,
111857cbf51SRichard Tran Mills   if each process wants only some of the values it should extract the ones it wants from the array.
112857cbf51SRichard Tran Mills 
1131cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetColumnSumsImaginaryPart()`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
114857cbf51SRichard Tran Mills @*/
MatGetColumnSumsRealPart(Mat A,PetscReal sums[])115d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnSumsRealPart(Mat A, PetscReal sums[])
116d71ae5a4SJacob Faibussowitsch {
117857cbf51SRichard Tran Mills   PetscFunctionBegin;
1189566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_SUM_REALPART, sums));
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120857cbf51SRichard Tran Mills }
121857cbf51SRichard Tran Mills 
122857cbf51SRichard Tran Mills /*@
123857cbf51SRichard Tran Mills   MatGetColumnSumsImaginaryPart - Gets the sums of the imaginary part of each column of a sparse or dense matrix.
124857cbf51SRichard Tran Mills 
125857cbf51SRichard Tran Mills   Input Parameter:
126857cbf51SRichard Tran Mills . A - the matrix
127857cbf51SRichard Tran Mills 
128857cbf51SRichard Tran Mills   Output Parameter:
129857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix
130857cbf51SRichard Tran Mills 
131857cbf51SRichard Tran Mills   Level: intermediate
132857cbf51SRichard Tran Mills 
13311a5261eSBarry Smith   Note:
134857cbf51SRichard Tran Mills   Each process has ALL the column sums after the call. Because of the way this is computed each process gets all the values,
135857cbf51SRichard Tran Mills   if each process wants only some of the values it should extract the ones it wants from the array.
136857cbf51SRichard Tran Mills 
1371cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetColumnSumsRealPart()`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
138857cbf51SRichard Tran Mills @*/
MatGetColumnSumsImaginaryPart(Mat A,PetscReal sums[])139d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnSumsImaginaryPart(Mat A, PetscReal sums[])
140d71ae5a4SJacob Faibussowitsch {
141857cbf51SRichard Tran Mills   PetscFunctionBegin;
1429566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_SUM_IMAGINARYPART, sums));
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144857cbf51SRichard Tran Mills }
145857cbf51SRichard Tran Mills 
146857cbf51SRichard Tran Mills /*@
147857cbf51SRichard Tran Mills   MatGetColumnSums - Gets the sums of each column of a sparse or dense matrix.
148857cbf51SRichard Tran Mills 
149857cbf51SRichard Tran Mills   Input Parameter:
150857cbf51SRichard Tran Mills . A - the matrix
151857cbf51SRichard Tran Mills 
152857cbf51SRichard Tran Mills   Output Parameter:
153857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix
154857cbf51SRichard Tran Mills 
155857cbf51SRichard Tran Mills   Level: intermediate
156857cbf51SRichard Tran Mills 
15711a5261eSBarry Smith   Note:
158857cbf51SRichard Tran Mills   Each process has ALL the column sums after the call. Because of the way this is computed each process gets all the values,
159857cbf51SRichard Tran Mills   if each process wants only some of the values it should extract the ones it wants from the array.
160857cbf51SRichard Tran Mills 
1611cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
162857cbf51SRichard Tran Mills @*/
MatGetColumnSums(Mat A,PetscScalar sums[])163d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnSums(Mat A, PetscScalar sums[])
164d71ae5a4SJacob Faibussowitsch {
165857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
166857cbf51SRichard Tran Mills   PetscInt   i, n;
167857cbf51SRichard Tran Mills   PetscReal *work;
168857cbf51SRichard Tran Mills #endif
169857cbf51SRichard Tran Mills 
170857cbf51SRichard Tran Mills   PetscFunctionBegin;
171857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
1729566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsRealPart(A, sums));
173857cbf51SRichard Tran Mills #else
1749566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &n));
1759566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(sums, n));
1769566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &work));
1779566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsRealPart(A, work));
178857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) sums[i] = work[i];
1799566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsImaginaryPart(A, work));
180857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) sums[i] += work[i] * PETSC_i;
1819566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
182857cbf51SRichard Tran Mills #endif
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
184857cbf51SRichard Tran Mills }
185857cbf51SRichard Tran Mills 
186857cbf51SRichard Tran Mills /*@
187857cbf51SRichard Tran Mills   MatGetColumnMeansRealPart - Gets the arithmetic means of the real part of each column of a sparse or dense matrix.
188857cbf51SRichard Tran Mills 
189857cbf51SRichard Tran Mills   Input Parameter:
190857cbf51SRichard Tran Mills . A - the matrix
191857cbf51SRichard Tran Mills 
192857cbf51SRichard Tran Mills   Output Parameter:
193fe59aa6dSJacob Faibussowitsch . means - an array as large as the TOTAL number of columns in the matrix
194857cbf51SRichard Tran Mills 
195857cbf51SRichard Tran Mills   Level: intermediate
196857cbf51SRichard Tran Mills 
19711a5261eSBarry Smith   Note:
198857cbf51SRichard Tran Mills   Each process has ALL the column means after the call. Because of the way this is computed each process gets all the values,
199857cbf51SRichard Tran Mills   if each process wants only some of the values it should extract the ones it wants from the array.
200857cbf51SRichard Tran Mills 
2011cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetColumnMeansImaginaryPart()`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
202857cbf51SRichard Tran Mills @*/
MatGetColumnMeansRealPart(Mat A,PetscReal means[])203d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnMeansRealPart(Mat A, PetscReal means[])
204d71ae5a4SJacob Faibussowitsch {
205857cbf51SRichard Tran Mills   PetscFunctionBegin;
2069566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_MEAN_REALPART, means));
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
208857cbf51SRichard Tran Mills }
209857cbf51SRichard Tran Mills 
210857cbf51SRichard Tran Mills /*@
211857cbf51SRichard Tran Mills   MatGetColumnMeansImaginaryPart - Gets the arithmetic means of the imaginary part of each column of a sparse or dense matrix.
212857cbf51SRichard Tran Mills 
213857cbf51SRichard Tran Mills   Input Parameter:
214857cbf51SRichard Tran Mills . A - the matrix
215857cbf51SRichard Tran Mills 
216857cbf51SRichard Tran Mills   Output Parameter:
217fe59aa6dSJacob Faibussowitsch . means - an array as large as the TOTAL number of columns in the matrix
218857cbf51SRichard Tran Mills 
219857cbf51SRichard Tran Mills   Level: intermediate
220857cbf51SRichard Tran Mills 
22111a5261eSBarry Smith   Note:
222857cbf51SRichard Tran Mills   Each process has ALL the column means after the call. Because of the way this is computed each process gets all the values,
223857cbf51SRichard Tran Mills   if each process wants only some of the values it should extract the ones it wants from the array.
224857cbf51SRichard Tran Mills 
2251cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetColumnMeansRealPart()`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
226857cbf51SRichard Tran Mills @*/
MatGetColumnMeansImaginaryPart(Mat A,PetscReal means[])227d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnMeansImaginaryPart(Mat A, PetscReal means[])
228d71ae5a4SJacob Faibussowitsch {
229857cbf51SRichard Tran Mills   PetscFunctionBegin;
2309566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_MEAN_IMAGINARYPART, means));
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
232857cbf51SRichard Tran Mills }
233857cbf51SRichard Tran Mills 
234857cbf51SRichard Tran Mills /*@
235857cbf51SRichard Tran Mills   MatGetColumnMeans - Gets the arithmetic means of each column of a sparse or dense matrix.
236857cbf51SRichard Tran Mills 
237857cbf51SRichard Tran Mills   Input Parameter:
238857cbf51SRichard Tran Mills . A - the matrix
239857cbf51SRichard Tran Mills 
240857cbf51SRichard Tran Mills   Output Parameter:
241857cbf51SRichard Tran Mills . means - an array as large as the TOTAL number of columns in the matrix
242857cbf51SRichard Tran Mills 
243857cbf51SRichard Tran Mills   Level: intermediate
244857cbf51SRichard Tran Mills 
24511a5261eSBarry Smith   Note:
246857cbf51SRichard Tran Mills   Each process has ALL the column means after the call. Because of the way this is computed each process gets all the values,
247857cbf51SRichard Tran Mills   if each process wants only some of the values it should extract the ones it wants from the array.
248857cbf51SRichard Tran Mills 
2491cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
250857cbf51SRichard Tran Mills @*/
MatGetColumnMeans(Mat A,PetscScalar means[])251d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnMeans(Mat A, PetscScalar means[])
252d71ae5a4SJacob Faibussowitsch {
253857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
254857cbf51SRichard Tran Mills   PetscInt   i, n;
255857cbf51SRichard Tran Mills   PetscReal *work;
256857cbf51SRichard Tran Mills #endif
257857cbf51SRichard Tran Mills 
258857cbf51SRichard Tran Mills   PetscFunctionBegin;
259857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
2609566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansRealPart(A, means));
261857cbf51SRichard Tran Mills #else
2629566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &n));
2639566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(means, n));
2649566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &work));
2659566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansRealPart(A, work));
266857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) means[i] = work[i];
2679566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansImaginaryPart(A, work));
268857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) means[i] += work[i] * PETSC_i;
2699566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
270857cbf51SRichard Tran Mills #endif
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272242f1d38SBarry Smith }
273242f1d38SBarry Smith 
274a873a8cdSSam Reynolds /*@
275a873a8cdSSam Reynolds   MatGetColumnReductions - Gets the reductions of each column of a sparse or dense matrix.
276a873a8cdSSam Reynolds 
27702024617SSatish Balay   Input Parameters:
278a873a8cdSSam Reynolds + A    - the matrix
27911a5261eSBarry Smith - type - A constant defined in `NormType` or `ReductionType`: `NORM_2`, `NORM_1`, `NORM_INFINITY`, `REDUCTION_SUM_REALPART`,
28011a5261eSBarry Smith           `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_MEAN_REALPART`, `REDUCTION_MEAN_IMAGINARYPART`
281a873a8cdSSam Reynolds 
282a873a8cdSSam Reynolds   Output Parameter:
283a873a8cdSSam Reynolds . reductions - an array as large as the TOTAL number of columns in the matrix
284a873a8cdSSam Reynolds 
285857cbf51SRichard Tran Mills   Level: developer
286a873a8cdSSam Reynolds 
28711a5261eSBarry Smith   Note:
288a873a8cdSSam Reynolds   Each process has ALL the column reductions after the call. Because of the way this is computed each process gets all the values,
289a873a8cdSSam Reynolds   if each process wants only some of the values it should extract the ones it wants from the array.
290a873a8cdSSam Reynolds 
291fe59aa6dSJacob Faibussowitsch   Developer Notes:
292857cbf51SRichard Tran Mills   This routine is primarily intended as a back-end.
29311a5261eSBarry Smith   `MatGetColumnNorms()`, `MatGetColumnSums()`, and `MatGetColumnMeans()` are implemented using this routine.
294a873a8cdSSam Reynolds 
2951cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `ReductionType`, `NormType`, `MatGetColumnNorms()`, `MatGetColumnSums()`, `MatGetColumnMeans()`
296a873a8cdSSam Reynolds @*/
MatGetColumnReductions(Mat A,PetscInt type,PetscReal reductions[])297d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnReductions(Mat A, PetscInt type, PetscReal reductions[])
298d71ae5a4SJacob Faibussowitsch {
299a873a8cdSSam Reynolds   PetscFunctionBegin;
300a873a8cdSSam Reynolds   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
301dbbe0bcdSBarry Smith   PetscUseTypeMethod(A, getcolumnreductions, type, reductions);
3023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
303a873a8cdSSam Reynolds }
304