xref: /petsc/src/mat/utils/getcolv.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
10925cdddSBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I   "petscmat.h"  I*/
30925cdddSBarry Smith 
40925cdddSBarry Smith /*@
582bf6240SBarry Smith    MatGetColumnVector - Gets the values from a given column of a matrix.
60925cdddSBarry Smith 
772257631SBarry Smith    Not Collective
8fee21e36SBarry Smith 
998a79cdbSBarry Smith    Input Parameters:
105ed6d96aSBarry Smith +  A - the matrix
115ed6d96aSBarry Smith .  yy - the vector
12b7e53d48SBarry Smith -  col - the column requested (in global numbering)
1398a79cdbSBarry Smith 
1415091d37SBarry Smith    Level: advanced
1515091d37SBarry Smith 
169d006df2SBarry Smith    Notes:
1727430b45SBarry Smith    If a `MatType` does not implement the operation, each processor for which this is called
1811a5261eSBarry Smith    gets the values for its rows using `MatGetRow()`.
199d006df2SBarry Smith 
203d81755aSBarry Smith    The vector must have the same parallel row layout as the matrix.
213d81755aSBarry Smith 
2282bf6240SBarry Smith    Contributed by: Denis Vanderstraeten
230925cdddSBarry Smith 
24*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetRow()`, `MatGetDiagonal()`, `MatMult()`
250925cdddSBarry Smith @*/
26d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnVector(Mat A, Vec yy, PetscInt col)
27d71ae5a4SJacob Faibussowitsch {
288aa348c1SBarry Smith   PetscScalar       *y;
29b3cc6726SBarry Smith   const PetscScalar *v;
3038baddfdSBarry Smith   PetscInt           i, j, nz, N, Rs, Re, rs, re;
3138baddfdSBarry Smith   const PetscInt    *idx;
320925cdddSBarry Smith 
330925cdddSBarry Smith   PetscFunctionBegin;
340700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
350700a824SBarry Smith   PetscValidHeaderSpecific(yy, VEC_CLASSID, 2);
369a971ab9SStefano Zampini   PetscValidLogicalCollectiveInt(A, col, 3);
3708401ef6SPierre Jolivet   PetscCheck(col >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Requested negative column: %" PetscInt_FMT, col);
389566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &N));
3908401ef6SPierre 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);
409566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &Rs, &Re));
419566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(yy, &rs, &re));
42aed4548fSBarry 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);
4382bf6240SBarry Smith 
44dbbe0bcdSBarry Smith   if (A->ops->getcolumnvector) PetscUseTypeMethod(A, getcolumnvector, yy, col);
45dbbe0bcdSBarry Smith   else {
469566063dSJacob Faibussowitsch     PetscCall(VecSet(yy, 0.0));
479566063dSJacob Faibussowitsch     PetscCall(VecGetArray(yy, &y));
489a971ab9SStefano Zampini     /* TODO for general matrices */
4982bf6240SBarry Smith     for (i = Rs; i < Re; i++) {
509566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A, i, &nz, &idx, &v));
5182bf6240SBarry Smith       if (nz && idx[0] <= col) {
5282bf6240SBarry Smith         /*
5382bf6240SBarry Smith           Should use faster search here
5482bf6240SBarry Smith         */
5582bf6240SBarry Smith         for (j = 0; j < nz; j++) {
5682bf6240SBarry Smith           if (idx[j] >= col) {
5782bf6240SBarry Smith             if (idx[j] == col) y[i - rs] = v[j];
5882bf6240SBarry Smith             break;
590925cdddSBarry Smith           }
600925cdddSBarry Smith         }
610925cdddSBarry Smith       }
629566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A, i, &nz, &idx, &v));
630925cdddSBarry Smith     }
649566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(yy, &y));
658d0534beSBarry Smith   }
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
670925cdddSBarry Smith }
68242f1d38SBarry Smith 
6986819fdcSBarry Smith /*@
700716a85fSBarry Smith    MatGetColumnNorms - Gets the norms of each column of a sparse or dense matrix.
71242f1d38SBarry Smith 
72d8d19677SJose E. Roman    Input Parameters:
73242f1d38SBarry Smith +  A - the matrix
7411a5261eSBarry Smith -  type - `NORM_2`, `NORM_1` or `NORM_INFINITY`
75242f1d38SBarry Smith 
76242f1d38SBarry Smith    Output Parameter:
77242f1d38SBarry Smith .  norms - an array as large as the TOTAL number of columns in the matrix
78242f1d38SBarry Smith 
79f6680f47SSatish Balay    Level: intermediate
80f6680f47SSatish Balay 
8111a5261eSBarry Smith    Note:
8295452b02SPatrick Sanan     Each process has ALL the column norms after the call. Because of the way this is computed each process gets all the values,
8386819fdcSBarry Smith     if each process wants only some of the values it should extract the ones it wants from the array.
8486819fdcSBarry Smith 
85*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `NormType`, `MatNorm()`
8686819fdcSBarry Smith @*/
87d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnNorms(Mat A, NormType type, PetscReal norms[])
88d71ae5a4SJacob Faibussowitsch {
89857cbf51SRichard Tran Mills   /* NOTE: MatGetColumnNorms() could simply be a macro that calls MatGetColumnReductions().
90857cbf51SRichard Tran Mills    * I've kept this as a function because it allows slightly more in the way of error checking,
91857cbf51SRichard Tran Mills    * erroring out if MatGetColumnNorms() is not called with a valid NormType. */
92242f1d38SBarry Smith 
93242f1d38SBarry Smith   PetscFunctionBegin;
94857cbf51SRichard Tran Mills   if (type == NORM_2 || type == NORM_1 || type == NORM_FROBENIUS || type == NORM_INFINITY || type == NORM_1_AND_2) {
959566063dSJacob Faibussowitsch     PetscCall(MatGetColumnReductions(A, type, norms));
96857cbf51SRichard Tran Mills   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Unknown NormType");
973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98a873a8cdSSam Reynolds }
99857cbf51SRichard Tran Mills 
100857cbf51SRichard Tran Mills /*@
101857cbf51SRichard Tran Mills    MatGetColumnSumsRealPart - Gets the sums of the real part of each column of a sparse or dense matrix.
102857cbf51SRichard Tran Mills 
103857cbf51SRichard Tran Mills    Input Parameter:
104857cbf51SRichard Tran Mills .  A - the matrix
105857cbf51SRichard Tran Mills 
106857cbf51SRichard Tran Mills    Output Parameter:
107857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
108857cbf51SRichard Tran Mills 
109857cbf51SRichard Tran Mills    Level: intermediate
110857cbf51SRichard Tran Mills 
11111a5261eSBarry Smith    Note:
112857cbf51SRichard 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,
113857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
114857cbf51SRichard Tran Mills 
115*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetColumnSumsImaginaryPart()`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
116857cbf51SRichard Tran Mills @*/
117d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnSumsRealPart(Mat A, PetscReal sums[])
118d71ae5a4SJacob Faibussowitsch {
119857cbf51SRichard Tran Mills   PetscFunctionBegin;
1209566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_SUM_REALPART, sums));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122857cbf51SRichard Tran Mills }
123857cbf51SRichard Tran Mills 
124857cbf51SRichard Tran Mills /*@
125857cbf51SRichard Tran Mills    MatGetColumnSumsImaginaryPart - Gets the sums of the imaginary part of each column of a sparse or dense matrix.
126857cbf51SRichard Tran Mills 
127857cbf51SRichard Tran Mills    Input Parameter:
128857cbf51SRichard Tran Mills .  A - the matrix
129857cbf51SRichard Tran Mills 
130857cbf51SRichard Tran Mills    Output Parameter:
131857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
132857cbf51SRichard Tran Mills 
133857cbf51SRichard Tran Mills    Level: intermediate
134857cbf51SRichard Tran Mills 
13511a5261eSBarry Smith    Note:
136857cbf51SRichard 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,
137857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
138857cbf51SRichard Tran Mills 
139*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetColumnSumsRealPart()`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
140857cbf51SRichard Tran Mills @*/
141d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnSumsImaginaryPart(Mat A, PetscReal sums[])
142d71ae5a4SJacob Faibussowitsch {
143857cbf51SRichard Tran Mills   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_SUM_IMAGINARYPART, sums));
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
146857cbf51SRichard Tran Mills }
147857cbf51SRichard Tran Mills 
148857cbf51SRichard Tran Mills /*@
149857cbf51SRichard Tran Mills    MatGetColumnSums - Gets the sums of each column of a sparse or dense matrix.
150857cbf51SRichard Tran Mills 
151857cbf51SRichard Tran Mills    Input Parameter:
152857cbf51SRichard Tran Mills .  A - the matrix
153857cbf51SRichard Tran Mills 
154857cbf51SRichard Tran Mills    Output Parameter:
155857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
156857cbf51SRichard Tran Mills 
157857cbf51SRichard Tran Mills    Level: intermediate
158857cbf51SRichard Tran Mills 
15911a5261eSBarry Smith    Note:
160857cbf51SRichard 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,
161857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
162857cbf51SRichard Tran Mills 
163*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
164857cbf51SRichard Tran Mills @*/
165d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnSums(Mat A, PetscScalar sums[])
166d71ae5a4SJacob Faibussowitsch {
167857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
168857cbf51SRichard Tran Mills   PetscInt   i, n;
169857cbf51SRichard Tran Mills   PetscReal *work;
170857cbf51SRichard Tran Mills #endif
171857cbf51SRichard Tran Mills 
172857cbf51SRichard Tran Mills   PetscFunctionBegin;
173857cbf51SRichard Tran Mills 
174857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
1759566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsRealPart(A, sums));
176857cbf51SRichard Tran Mills #else
1779566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &n));
1789566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(sums, n));
1799566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &work));
1809566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsRealPart(A, work));
181857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) sums[i] = work[i];
1829566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsImaginaryPart(A, work));
183857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) sums[i] += work[i] * PETSC_i;
1849566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
185857cbf51SRichard Tran Mills #endif
1863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
187857cbf51SRichard Tran Mills }
188857cbf51SRichard Tran Mills 
189857cbf51SRichard Tran Mills /*@
190857cbf51SRichard Tran Mills    MatGetColumnMeansRealPart - Gets the arithmetic means of the real part of each column of a sparse or dense matrix.
191857cbf51SRichard Tran Mills 
192857cbf51SRichard Tran Mills    Input Parameter:
193857cbf51SRichard Tran Mills .  A - the matrix
194857cbf51SRichard Tran Mills 
195857cbf51SRichard Tran Mills    Output Parameter:
196857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
197857cbf51SRichard Tran Mills 
198857cbf51SRichard Tran Mills    Level: intermediate
199857cbf51SRichard Tran Mills 
20011a5261eSBarry Smith    Note:
201857cbf51SRichard 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,
202857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
203857cbf51SRichard Tran Mills 
204*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetColumnMeansImaginaryPart()`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
205857cbf51SRichard Tran Mills @*/
206d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnMeansRealPart(Mat A, PetscReal means[])
207d71ae5a4SJacob Faibussowitsch {
208857cbf51SRichard Tran Mills   PetscFunctionBegin;
2099566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_MEAN_REALPART, means));
2103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
211857cbf51SRichard Tran Mills }
212857cbf51SRichard Tran Mills 
213857cbf51SRichard Tran Mills /*@
214857cbf51SRichard Tran Mills    MatGetColumnMeansImaginaryPart - Gets the arithmetic means of the imaginary part of each column of a sparse or dense matrix.
215857cbf51SRichard Tran Mills 
216857cbf51SRichard Tran Mills    Input Parameter:
217857cbf51SRichard Tran Mills .  A - the matrix
218857cbf51SRichard Tran Mills 
219857cbf51SRichard Tran Mills    Output Parameter:
220857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
221857cbf51SRichard Tran Mills 
222857cbf51SRichard Tran Mills    Level: intermediate
223857cbf51SRichard Tran Mills 
22411a5261eSBarry Smith    Note:
225857cbf51SRichard 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,
226857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
227857cbf51SRichard Tran Mills 
228*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatGetColumnMeansRealPart()`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
229857cbf51SRichard Tran Mills @*/
230d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnMeansImaginaryPart(Mat A, PetscReal means[])
231d71ae5a4SJacob Faibussowitsch {
232857cbf51SRichard Tran Mills   PetscFunctionBegin;
2339566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_MEAN_IMAGINARYPART, means));
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
235857cbf51SRichard Tran Mills }
236857cbf51SRichard Tran Mills 
237857cbf51SRichard Tran Mills /*@
238857cbf51SRichard Tran Mills    MatGetColumnMeans - Gets the arithmetic means of each column of a sparse or dense matrix.
239857cbf51SRichard Tran Mills 
240857cbf51SRichard Tran Mills    Input Parameter:
241857cbf51SRichard Tran Mills .  A - the matrix
242857cbf51SRichard Tran Mills 
243857cbf51SRichard Tran Mills    Output Parameter:
244857cbf51SRichard Tran Mills .  means - an array as large as the TOTAL number of columns in the matrix
245857cbf51SRichard Tran Mills 
246857cbf51SRichard Tran Mills    Level: intermediate
247857cbf51SRichard Tran Mills 
24811a5261eSBarry Smith    Note:
249857cbf51SRichard 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,
250857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
251857cbf51SRichard Tran Mills 
252*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
253857cbf51SRichard Tran Mills @*/
254d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnMeans(Mat A, PetscScalar means[])
255d71ae5a4SJacob Faibussowitsch {
256857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
257857cbf51SRichard Tran Mills   PetscInt   i, n;
258857cbf51SRichard Tran Mills   PetscReal *work;
259857cbf51SRichard Tran Mills #endif
260857cbf51SRichard Tran Mills 
261857cbf51SRichard Tran Mills   PetscFunctionBegin;
262857cbf51SRichard Tran Mills 
263857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
2649566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansRealPart(A, means));
265857cbf51SRichard Tran Mills #else
2669566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &n));
2679566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(means, n));
2689566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &work));
2699566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansRealPart(A, work));
270857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) means[i] = work[i];
2719566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansImaginaryPart(A, work));
272857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) means[i] += work[i] * PETSC_i;
2739566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
274857cbf51SRichard Tran Mills #endif
2753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
276242f1d38SBarry Smith }
277242f1d38SBarry Smith 
278a873a8cdSSam Reynolds /*@
279a873a8cdSSam Reynolds     MatGetColumnReductions - Gets the reductions of each column of a sparse or dense matrix.
280a873a8cdSSam Reynolds 
28102024617SSatish Balay   Input Parameters:
282a873a8cdSSam Reynolds +  A - the matrix
28311a5261eSBarry Smith -  type - A constant defined in `NormType` or `ReductionType`: `NORM_2`, `NORM_1`, `NORM_INFINITY`, `REDUCTION_SUM_REALPART`,
28411a5261eSBarry Smith           `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_MEAN_REALPART`, `REDUCTION_MEAN_IMAGINARYPART`
285a873a8cdSSam Reynolds 
286a873a8cdSSam Reynolds   Output Parameter:
287a873a8cdSSam Reynolds .  reductions - an array as large as the TOTAL number of columns in the matrix
288a873a8cdSSam Reynolds 
289857cbf51SRichard Tran Mills    Level: developer
290a873a8cdSSam Reynolds 
29111a5261eSBarry Smith    Note:
292a873a8cdSSam Reynolds     Each process has ALL the column reductions after the call. Because of the way this is computed each process gets all the values,
293a873a8cdSSam Reynolds     if each process wants only some of the values it should extract the ones it wants from the array.
294a873a8cdSSam Reynolds 
295a873a8cdSSam Reynolds   Developer Note:
296857cbf51SRichard Tran Mills     This routine is primarily intended as a back-end.
29711a5261eSBarry Smith     `MatGetColumnNorms()`, `MatGetColumnSums()`, and `MatGetColumnMeans()` are implemented using this routine.
298a873a8cdSSam Reynolds 
299*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `ReductionType`, `NormType`, `MatGetColumnNorms()`, `MatGetColumnSums()`, `MatGetColumnMeans()`
300a873a8cdSSam Reynolds @*/
301d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnReductions(Mat A, PetscInt type, PetscReal reductions[])
302d71ae5a4SJacob Faibussowitsch {
303a873a8cdSSam Reynolds   PetscFunctionBegin;
304a873a8cdSSam Reynolds   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
305dbbe0bcdSBarry Smith   PetscUseTypeMethod(A, getcolumnreductions, type, reductions);
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307a873a8cdSSam Reynolds }
308