xref: /petsc/src/mat/utils/getcolv.c (revision 11a5261e40035b7c793f2783a2ba6c7cd4f3b077)
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:
179a971ab9SStefano Zampini    If a Mat type does not implement the operation, each processor for which this is called
18*11a5261eSBarry 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 
24db781477SPatrick Sanan .seealso: `MatGetRow()`, `MatGetDiagonal()`, `MatMult()`
250925cdddSBarry Smith @*/
269371c9d4SSatish Balay PetscErrorCode MatGetColumnVector(Mat A, Vec yy, PetscInt col) {
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   }
650925cdddSBarry Smith   PetscFunctionReturn(0);
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
73*11a5261eSBarry 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 
80*11a5261eSBarry 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 
84db781477SPatrick Sanan .seealso: `NormType`, `MatNorm()`
8586819fdcSBarry Smith @*/
869371c9d4SSatish Balay PetscErrorCode MatGetColumnNorms(Mat A, NormType type, PetscReal norms[]) {
87857cbf51SRichard Tran Mills   /* NOTE: MatGetColumnNorms() could simply be a macro that calls MatGetColumnReductions().
88857cbf51SRichard Tran Mills    * I've kept this as a function because it allows slightly more in the way of error checking,
89857cbf51SRichard Tran Mills    * erroring out if MatGetColumnNorms() is not called with a valid NormType. */
90242f1d38SBarry Smith 
91242f1d38SBarry Smith   PetscFunctionBegin;
92857cbf51SRichard Tran Mills   if (type == NORM_2 || type == NORM_1 || type == NORM_FROBENIUS || type == NORM_INFINITY || type == NORM_1_AND_2) {
939566063dSJacob Faibussowitsch     PetscCall(MatGetColumnReductions(A, type, norms));
94857cbf51SRichard Tran Mills   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Unknown NormType");
95857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
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 
109*11a5261eSBarry 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 
113db781477SPatrick Sanan .seealso: `MatGetColumnSumsImaginaryPart()`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
114857cbf51SRichard Tran Mills @*/
1159371c9d4SSatish Balay PetscErrorCode MatGetColumnSumsRealPart(Mat A, PetscReal sums[]) {
116857cbf51SRichard Tran Mills   PetscFunctionBegin;
1179566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_SUM_REALPART, sums));
118857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
119857cbf51SRichard Tran Mills }
120857cbf51SRichard Tran Mills 
121857cbf51SRichard Tran Mills /*@
122857cbf51SRichard Tran Mills    MatGetColumnSumsImaginaryPart - Gets the sums of the imaginary part of each column of a sparse or dense matrix.
123857cbf51SRichard Tran Mills 
124857cbf51SRichard Tran Mills    Input Parameter:
125857cbf51SRichard Tran Mills .  A - the matrix
126857cbf51SRichard Tran Mills 
127857cbf51SRichard Tran Mills    Output Parameter:
128857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
129857cbf51SRichard Tran Mills 
130857cbf51SRichard Tran Mills    Level: intermediate
131857cbf51SRichard Tran Mills 
132*11a5261eSBarry Smith    Note:
133857cbf51SRichard 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,
134857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
135857cbf51SRichard Tran Mills 
136db781477SPatrick Sanan .seealso: `MatGetColumnSumsRealPart()`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
137857cbf51SRichard Tran Mills @*/
1389371c9d4SSatish Balay PetscErrorCode MatGetColumnSumsImaginaryPart(Mat A, PetscReal sums[]) {
139857cbf51SRichard Tran Mills   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_SUM_IMAGINARYPART, sums));
141857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
142857cbf51SRichard Tran Mills }
143857cbf51SRichard Tran Mills 
144857cbf51SRichard Tran Mills /*@
145857cbf51SRichard Tran Mills    MatGetColumnSums - Gets the sums of each column of a sparse or dense matrix.
146857cbf51SRichard Tran Mills 
147857cbf51SRichard Tran Mills    Input Parameter:
148857cbf51SRichard Tran Mills .  A - the matrix
149857cbf51SRichard Tran Mills 
150857cbf51SRichard Tran Mills    Output Parameter:
151857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
152857cbf51SRichard Tran Mills 
153857cbf51SRichard Tran Mills    Level: intermediate
154857cbf51SRichard Tran Mills 
155*11a5261eSBarry Smith    Note:
156857cbf51SRichard 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,
157857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
158857cbf51SRichard Tran Mills 
159db781477SPatrick Sanan .seealso: `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
160857cbf51SRichard Tran Mills @*/
1619371c9d4SSatish Balay PetscErrorCode MatGetColumnSums(Mat A, PetscScalar sums[]) {
162857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
163857cbf51SRichard Tran Mills   PetscInt   i, n;
164857cbf51SRichard Tran Mills   PetscReal *work;
165857cbf51SRichard Tran Mills #endif
166857cbf51SRichard Tran Mills 
167857cbf51SRichard Tran Mills   PetscFunctionBegin;
168857cbf51SRichard Tran Mills 
169857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
1709566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsRealPart(A, sums));
171857cbf51SRichard Tran Mills #else
1729566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &n));
1739566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(sums, n));
1749566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &work));
1759566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsRealPart(A, work));
176857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) sums[i] = work[i];
1779566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsImaginaryPart(A, work));
178857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) sums[i] += work[i] * PETSC_i;
1799566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
180857cbf51SRichard Tran Mills #endif
181857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
182857cbf51SRichard Tran Mills }
183857cbf51SRichard Tran Mills 
184857cbf51SRichard Tran Mills /*@
185857cbf51SRichard Tran Mills    MatGetColumnMeansRealPart - Gets the arithmetic means of the real part of each column of a sparse or dense matrix.
186857cbf51SRichard Tran Mills 
187857cbf51SRichard Tran Mills    Input Parameter:
188857cbf51SRichard Tran Mills .  A - the matrix
189857cbf51SRichard Tran Mills 
190857cbf51SRichard Tran Mills    Output Parameter:
191857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
192857cbf51SRichard Tran Mills 
193857cbf51SRichard Tran Mills    Level: intermediate
194857cbf51SRichard Tran Mills 
195*11a5261eSBarry Smith    Note:
196857cbf51SRichard 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,
197857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
198857cbf51SRichard Tran Mills 
199db781477SPatrick Sanan .seealso: `MatGetColumnMeansImaginaryPart()`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
200857cbf51SRichard Tran Mills @*/
2019371c9d4SSatish Balay PetscErrorCode MatGetColumnMeansRealPart(Mat A, PetscReal means[]) {
202857cbf51SRichard Tran Mills   PetscFunctionBegin;
2039566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_MEAN_REALPART, means));
204857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
205857cbf51SRichard Tran Mills }
206857cbf51SRichard Tran Mills 
207857cbf51SRichard Tran Mills /*@
208857cbf51SRichard Tran Mills    MatGetColumnMeansImaginaryPart - Gets the arithmetic means of the imaginary part of each column of a sparse or dense matrix.
209857cbf51SRichard Tran Mills 
210857cbf51SRichard Tran Mills    Input Parameter:
211857cbf51SRichard Tran Mills .  A - the matrix
212857cbf51SRichard Tran Mills 
213857cbf51SRichard Tran Mills    Output Parameter:
214857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
215857cbf51SRichard Tran Mills 
216857cbf51SRichard Tran Mills    Level: intermediate
217857cbf51SRichard Tran Mills 
218*11a5261eSBarry Smith    Note:
219857cbf51SRichard 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,
220857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
221857cbf51SRichard Tran Mills 
222db781477SPatrick Sanan .seealso: `MatGetColumnMeansRealPart()`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
223857cbf51SRichard Tran Mills @*/
2249371c9d4SSatish Balay PetscErrorCode MatGetColumnMeansImaginaryPart(Mat A, PetscReal means[]) {
225857cbf51SRichard Tran Mills   PetscFunctionBegin;
2269566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A, REDUCTION_MEAN_IMAGINARYPART, means));
227857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
228857cbf51SRichard Tran Mills }
229857cbf51SRichard Tran Mills 
230857cbf51SRichard Tran Mills /*@
231857cbf51SRichard Tran Mills    MatGetColumnMeans - Gets the arithmetic means of each column of a sparse or dense matrix.
232857cbf51SRichard Tran Mills 
233857cbf51SRichard Tran Mills    Input Parameter:
234857cbf51SRichard Tran Mills .  A - the matrix
235857cbf51SRichard Tran Mills 
236857cbf51SRichard Tran Mills    Output Parameter:
237857cbf51SRichard Tran Mills .  means - an array as large as the TOTAL number of columns in the matrix
238857cbf51SRichard Tran Mills 
239857cbf51SRichard Tran Mills    Level: intermediate
240857cbf51SRichard Tran Mills 
241*11a5261eSBarry Smith    Note:
242857cbf51SRichard 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,
243857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
244857cbf51SRichard Tran Mills 
245db781477SPatrick Sanan .seealso: `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()`
246857cbf51SRichard Tran Mills @*/
2479371c9d4SSatish Balay PetscErrorCode MatGetColumnMeans(Mat A, PetscScalar means[]) {
248857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
249857cbf51SRichard Tran Mills   PetscInt   i, n;
250857cbf51SRichard Tran Mills   PetscReal *work;
251857cbf51SRichard Tran Mills #endif
252857cbf51SRichard Tran Mills 
253857cbf51SRichard Tran Mills   PetscFunctionBegin;
254857cbf51SRichard Tran Mills 
255857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
2569566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansRealPart(A, means));
257857cbf51SRichard Tran Mills #else
2589566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &n));
2599566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(means, n));
2609566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n, &work));
2619566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansRealPart(A, work));
262857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) means[i] = work[i];
2639566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansImaginaryPart(A, work));
264857cbf51SRichard Tran Mills   for (i = 0; i < n; i++) means[i] += work[i] * PETSC_i;
2659566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
266857cbf51SRichard Tran Mills #endif
267242f1d38SBarry Smith   PetscFunctionReturn(0);
268242f1d38SBarry Smith }
269242f1d38SBarry Smith 
270a873a8cdSSam Reynolds /*@
271a873a8cdSSam Reynolds     MatGetColumnReductions - Gets the reductions of each column of a sparse or dense matrix.
272a873a8cdSSam Reynolds 
27302024617SSatish Balay   Input Parameters:
274a873a8cdSSam Reynolds +  A - the matrix
275*11a5261eSBarry Smith -  type - A constant defined in `NormType` or `ReductionType`: `NORM_2`, `NORM_1`, `NORM_INFINITY`, `REDUCTION_SUM_REALPART`,
276*11a5261eSBarry Smith           `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_MEAN_REALPART`, `REDUCTION_MEAN_IMAGINARYPART`
277a873a8cdSSam Reynolds 
278a873a8cdSSam Reynolds   Output Parameter:
279a873a8cdSSam Reynolds .  reductions - an array as large as the TOTAL number of columns in the matrix
280a873a8cdSSam Reynolds 
281857cbf51SRichard Tran Mills    Level: developer
282a873a8cdSSam Reynolds 
283*11a5261eSBarry Smith    Note:
284a873a8cdSSam Reynolds     Each process has ALL the column reductions after the call. Because of the way this is computed each process gets all the values,
285a873a8cdSSam Reynolds     if each process wants only some of the values it should extract the ones it wants from the array.
286a873a8cdSSam Reynolds 
287a873a8cdSSam Reynolds   Developer Note:
288857cbf51SRichard Tran Mills     This routine is primarily intended as a back-end.
289*11a5261eSBarry Smith     `MatGetColumnNorms()`, `MatGetColumnSums()`, and `MatGetColumnMeans()` are implemented using this routine.
290a873a8cdSSam Reynolds 
291db781477SPatrick Sanan .seealso: `ReductionType`, `NormType`, `MatGetColumnNorms()`, `MatGetColumnSums()`, `MatGetColumnMeans()`
292a873a8cdSSam Reynolds @*/
2939371c9d4SSatish Balay PetscErrorCode MatGetColumnReductions(Mat A, PetscInt type, PetscReal reductions[]) {
294a873a8cdSSam Reynolds   PetscFunctionBegin;
295a873a8cdSSam Reynolds   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
296dbbe0bcdSBarry Smith   PetscUseTypeMethod(A, getcolumnreductions, type, reductions);
297a873a8cdSSam Reynolds   PetscFunctionReturn(0);
298a873a8cdSSam Reynolds }
299