xref: /petsc/src/mat/utils/getcolv.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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
189a971ab9SStefano Zampini    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 
249a971ab9SStefano Zampini .seealso: MatGetRow(), MatGetDiagonal(), MatMult()
2515091d37SBarry Smith 
260925cdddSBarry Smith @*/
277087cfbeSBarry Smith PetscErrorCode  MatGetColumnVector(Mat A,Vec yy,PetscInt col)
280925cdddSBarry Smith {
298aa348c1SBarry Smith   PetscScalar       *y;
30b3cc6726SBarry Smith   const PetscScalar *v;
3138baddfdSBarry Smith   PetscInt          i,j,nz,N,Rs,Re,rs,re;
3238baddfdSBarry Smith   const PetscInt    *idx;
330925cdddSBarry Smith 
340925cdddSBarry Smith   PetscFunctionBegin;
350700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
360700a824SBarry Smith   PetscValidHeaderSpecific(yy,VEC_CLASSID,2);
379a971ab9SStefano Zampini   PetscValidLogicalCollectiveInt(A,col,3);
38*08401ef6SPierre Jolivet   PetscCheck(col >= 0,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Requested negative column: %" PetscInt_FMT,col);
399566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,NULL,&N));
40*08401ef6SPierre 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);
419566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&Rs,&Re));
429566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(yy,&rs,&re));
432c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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);
4482bf6240SBarry Smith 
458d0534beSBarry Smith   if (A->ops->getcolumnvector) {
469566063dSJacob Faibussowitsch     PetscCall((*A->ops->getcolumnvector)(A,yy,col));
478d0534beSBarry Smith   } else {
489566063dSJacob Faibussowitsch     PetscCall(VecSet(yy,0.0));
499566063dSJacob Faibussowitsch     PetscCall(VecGetArray(yy,&y));
509a971ab9SStefano Zampini     /* TODO for general matrices */
5182bf6240SBarry Smith     for (i=Rs; i<Re; i++) {
529566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A,i,&nz,&idx,&v));
5382bf6240SBarry Smith       if (nz && idx[0] <= col) {
5482bf6240SBarry Smith         /*
5582bf6240SBarry Smith           Should use faster search here
5682bf6240SBarry Smith         */
5782bf6240SBarry Smith         for (j=0; j<nz; j++) {
5882bf6240SBarry Smith           if (idx[j] >= col) {
5982bf6240SBarry Smith             if (idx[j] == col) y[i-rs] = v[j];
6082bf6240SBarry Smith             break;
610925cdddSBarry Smith           }
620925cdddSBarry Smith         }
630925cdddSBarry Smith       }
649566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A,i,&nz,&idx,&v));
650925cdddSBarry Smith     }
669566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(yy,&y));
678d0534beSBarry Smith   }
680925cdddSBarry Smith   PetscFunctionReturn(0);
690925cdddSBarry Smith }
70242f1d38SBarry Smith 
7186819fdcSBarry Smith /*@
720716a85fSBarry Smith    MatGetColumnNorms - Gets the norms of each column of a sparse or dense matrix.
73242f1d38SBarry Smith 
74d8d19677SJose E. Roman    Input Parameters:
75242f1d38SBarry Smith +  A - the matrix
76242f1d38SBarry Smith -  type - NORM_2, NORM_1 or NORM_INFINITY
77242f1d38SBarry Smith 
78242f1d38SBarry Smith    Output Parameter:
79242f1d38SBarry Smith .  norms - an array as large as the TOTAL number of columns in the matrix
80242f1d38SBarry Smith 
81f6680f47SSatish Balay    Level: intermediate
82f6680f47SSatish Balay 
8395452b02SPatrick Sanan    Notes:
8495452b02SPatrick Sanan     Each process has ALL the column norms after the call. Because of the way this is computed each process gets all the values,
8586819fdcSBarry Smith     if each process wants only some of the values it should extract the ones it wants from the array.
8686819fdcSBarry Smith 
8740b1048aSBarry Smith .seealso: NormType, MatNorm()
8886819fdcSBarry Smith 
8986819fdcSBarry Smith @*/
906b9dab40SJed Brown PetscErrorCode MatGetColumnNorms(Mat A,NormType type,PetscReal norms[])
91242f1d38SBarry Smith {
92857cbf51SRichard Tran Mills   /* NOTE: MatGetColumnNorms() could simply be a macro that calls MatGetColumnReductions().
93857cbf51SRichard Tran Mills    * I've kept this as a function because it allows slightly more in the way of error checking,
94857cbf51SRichard Tran Mills    * erroring out if MatGetColumnNorms() is not called with a valid NormType. */
95242f1d38SBarry Smith 
96242f1d38SBarry Smith   PetscFunctionBegin;
97857cbf51SRichard Tran Mills   if (type == NORM_2 || type == NORM_1 || type == NORM_FROBENIUS || type == NORM_INFINITY || type == NORM_1_AND_2) {
989566063dSJacob Faibussowitsch     PetscCall(MatGetColumnReductions(A,type,norms));
99857cbf51SRichard Tran Mills   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType");
100857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
101a873a8cdSSam Reynolds }
102857cbf51SRichard Tran Mills 
103857cbf51SRichard Tran Mills /*@
104857cbf51SRichard Tran Mills    MatGetColumnSumsRealPart - Gets the sums of the real part of each column of a sparse or dense matrix.
105857cbf51SRichard Tran Mills 
106857cbf51SRichard Tran Mills    Input Parameter:
107857cbf51SRichard Tran Mills .  A - the matrix
108857cbf51SRichard Tran Mills 
109857cbf51SRichard Tran Mills    Output Parameter:
110857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
111857cbf51SRichard Tran Mills 
112857cbf51SRichard Tran Mills    Level: intermediate
113857cbf51SRichard Tran Mills 
114857cbf51SRichard Tran Mills    Notes:
115857cbf51SRichard 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,
116857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
117857cbf51SRichard Tran Mills 
118857cbf51SRichard Tran Mills .seealso: MatGetColumnSumsImaginaryPart(), VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions()
119857cbf51SRichard Tran Mills 
120857cbf51SRichard Tran Mills @*/
121857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSumsRealPart(Mat A,PetscReal sums[])
122857cbf51SRichard Tran Mills {
123857cbf51SRichard Tran Mills   PetscFunctionBegin;
1249566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A,REDUCTION_SUM_REALPART,sums));
125857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
126857cbf51SRichard Tran Mills }
127857cbf51SRichard Tran Mills 
128857cbf51SRichard Tran Mills /*@
129857cbf51SRichard Tran Mills    MatGetColumnSumsImaginaryPart - Gets the sums of the imaginary part of each column of a sparse or dense matrix.
130857cbf51SRichard Tran Mills 
131857cbf51SRichard Tran Mills    Input Parameter:
132857cbf51SRichard Tran Mills .  A - the matrix
133857cbf51SRichard Tran Mills 
134857cbf51SRichard Tran Mills    Output Parameter:
135857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
136857cbf51SRichard Tran Mills 
137857cbf51SRichard Tran Mills    Level: intermediate
138857cbf51SRichard Tran Mills 
139857cbf51SRichard Tran Mills    Notes:
140857cbf51SRichard 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,
141857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
142857cbf51SRichard Tran Mills 
143857cbf51SRichard Tran Mills .seealso: MatGetColumnSumsRealPart(), VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions()
144857cbf51SRichard Tran Mills 
145857cbf51SRichard Tran Mills @*/
146857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSumsImaginaryPart(Mat A,PetscReal sums[])
147857cbf51SRichard Tran Mills {
148857cbf51SRichard Tran Mills   PetscFunctionBegin;
1499566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A,REDUCTION_SUM_IMAGINARYPART,sums));
150857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
151857cbf51SRichard Tran Mills }
152857cbf51SRichard Tran Mills 
153857cbf51SRichard Tran Mills /*@
154857cbf51SRichard Tran Mills    MatGetColumnSums - Gets the sums of each column of a sparse or dense matrix.
155857cbf51SRichard Tran Mills 
156857cbf51SRichard Tran Mills    Input Parameter:
157857cbf51SRichard Tran Mills .  A - the matrix
158857cbf51SRichard Tran Mills 
159857cbf51SRichard Tran Mills    Output Parameter:
160857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
161857cbf51SRichard Tran Mills 
162857cbf51SRichard Tran Mills    Level: intermediate
163857cbf51SRichard Tran Mills 
164857cbf51SRichard Tran Mills    Notes:
165857cbf51SRichard 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,
166857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
167857cbf51SRichard Tran Mills 
168857cbf51SRichard Tran Mills .seealso: VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions()
169857cbf51SRichard Tran Mills 
170857cbf51SRichard Tran Mills @*/
171857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSums(Mat A,PetscScalar sums[])
172857cbf51SRichard Tran Mills {
173857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
174857cbf51SRichard Tran Mills   PetscInt       i,n;
175857cbf51SRichard Tran Mills   PetscReal      *work;
176857cbf51SRichard Tran Mills #endif
177857cbf51SRichard Tran Mills 
178857cbf51SRichard Tran Mills   PetscFunctionBegin;
179857cbf51SRichard Tran Mills 
180857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
1819566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsRealPart(A,sums));
182857cbf51SRichard Tran Mills #else
1839566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,NULL,&n));
1849566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(sums,n));
1859566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n,&work));
1869566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsRealPart(A,work));
187857cbf51SRichard Tran Mills   for (i=0; i<n; i++) sums[i] = work[i];
1889566063dSJacob Faibussowitsch   PetscCall(MatGetColumnSumsImaginaryPart(A,work));
189857cbf51SRichard Tran Mills   for (i=0; i<n; i++) sums[i] += work[i]*PETSC_i;
1909566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
191857cbf51SRichard Tran Mills #endif
192857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
193857cbf51SRichard Tran Mills }
194857cbf51SRichard Tran Mills 
195857cbf51SRichard Tran Mills /*@
196857cbf51SRichard Tran Mills    MatGetColumnMeansRealPart - Gets the arithmetic means of the real part of each column of a sparse or dense matrix.
197857cbf51SRichard Tran Mills 
198857cbf51SRichard Tran Mills    Input Parameter:
199857cbf51SRichard Tran Mills .  A - the matrix
200857cbf51SRichard Tran Mills 
201857cbf51SRichard Tran Mills    Output Parameter:
202857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
203857cbf51SRichard Tran Mills 
204857cbf51SRichard Tran Mills    Level: intermediate
205857cbf51SRichard Tran Mills 
206857cbf51SRichard Tran Mills    Notes:
207857cbf51SRichard 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,
208857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
209857cbf51SRichard Tran Mills 
210857cbf51SRichard Tran Mills .seealso: MatGetColumnMeansImaginaryPart(), VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions()
211857cbf51SRichard Tran Mills 
212857cbf51SRichard Tran Mills @*/
213857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeansRealPart(Mat A,PetscReal means[])
214857cbf51SRichard Tran Mills {
215857cbf51SRichard Tran Mills   PetscFunctionBegin;
2169566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A,REDUCTION_MEAN_REALPART,means));
217857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
218857cbf51SRichard Tran Mills }
219857cbf51SRichard Tran Mills 
220857cbf51SRichard Tran Mills /*@
221857cbf51SRichard Tran Mills    MatGetColumnMeansImaginaryPart - Gets the arithmetic means of the imaginary part of each column of a sparse or dense matrix.
222857cbf51SRichard Tran Mills 
223857cbf51SRichard Tran Mills    Input Parameter:
224857cbf51SRichard Tran Mills .  A - the matrix
225857cbf51SRichard Tran Mills 
226857cbf51SRichard Tran Mills    Output Parameter:
227857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
228857cbf51SRichard Tran Mills 
229857cbf51SRichard Tran Mills    Level: intermediate
230857cbf51SRichard Tran Mills 
231857cbf51SRichard Tran Mills    Notes:
232857cbf51SRichard 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,
233857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
234857cbf51SRichard Tran Mills 
235857cbf51SRichard Tran Mills .seealso: MatGetColumnMeansRealPart(), VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions()
236857cbf51SRichard Tran Mills 
237857cbf51SRichard Tran Mills @*/
238857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeansImaginaryPart(Mat A,PetscReal means[])
239857cbf51SRichard Tran Mills {
240857cbf51SRichard Tran Mills   PetscFunctionBegin;
2419566063dSJacob Faibussowitsch   PetscCall(MatGetColumnReductions(A,REDUCTION_MEAN_IMAGINARYPART,means));
242857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
243857cbf51SRichard Tran Mills }
244857cbf51SRichard Tran Mills 
245857cbf51SRichard Tran Mills /*@
246857cbf51SRichard Tran Mills    MatGetColumnMeans - Gets the arithmetic means of each column of a sparse or dense matrix.
247857cbf51SRichard Tran Mills 
248857cbf51SRichard Tran Mills    Input Parameter:
249857cbf51SRichard Tran Mills .  A - the matrix
250857cbf51SRichard Tran Mills 
251857cbf51SRichard Tran Mills    Output Parameter:
252857cbf51SRichard Tran Mills .  means - an array as large as the TOTAL number of columns in the matrix
253857cbf51SRichard Tran Mills 
254857cbf51SRichard Tran Mills    Level: intermediate
255857cbf51SRichard Tran Mills 
256857cbf51SRichard Tran Mills    Notes:
257857cbf51SRichard 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,
258857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
259857cbf51SRichard Tran Mills 
260857cbf51SRichard Tran Mills .seealso: VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions()
261857cbf51SRichard Tran Mills 
262857cbf51SRichard Tran Mills @*/
263857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeans(Mat A,PetscScalar means[])
264857cbf51SRichard Tran Mills {
265857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
266857cbf51SRichard Tran Mills   PetscInt       i,n;
267857cbf51SRichard Tran Mills   PetscReal      *work;
268857cbf51SRichard Tran Mills #endif
269857cbf51SRichard Tran Mills 
270857cbf51SRichard Tran Mills   PetscFunctionBegin;
271857cbf51SRichard Tran Mills 
272857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
2739566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansRealPart(A,means));
274857cbf51SRichard Tran Mills #else
2759566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,NULL,&n));
2769566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(means,n));
2779566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(n,&work));
2789566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansRealPart(A,work));
279857cbf51SRichard Tran Mills   for (i=0; i<n; i++) means[i] = work[i];
2809566063dSJacob Faibussowitsch   PetscCall(MatGetColumnMeansImaginaryPart(A,work));
281857cbf51SRichard Tran Mills   for (i=0; i<n; i++) means[i] += work[i]*PETSC_i;
2829566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
283857cbf51SRichard Tran Mills #endif
284242f1d38SBarry Smith   PetscFunctionReturn(0);
285242f1d38SBarry Smith }
286242f1d38SBarry Smith 
287a873a8cdSSam Reynolds /*@
288a873a8cdSSam Reynolds     MatGetColumnReductions - Gets the reductions of each column of a sparse or dense matrix.
289a873a8cdSSam Reynolds 
29002024617SSatish Balay   Input Parameters:
291a873a8cdSSam Reynolds +  A - the matrix
292857cbf51SRichard Tran Mills -  type - A constant defined in NormType or ReductionType: NORM_2, NORM_1, NORM_INFINITY, REDUCTION_SUM_REALPART,
293857cbf51SRichard Tran Mills           REDUCTION_SUM_IMAGINARYPART, REDUCTION_MEAN_REALPART, REDUCTION_MEAN_IMAGINARYPART
294a873a8cdSSam Reynolds 
295a873a8cdSSam Reynolds   Output Parameter:
296a873a8cdSSam Reynolds .  reductions - an array as large as the TOTAL number of columns in the matrix
297a873a8cdSSam Reynolds 
298857cbf51SRichard Tran Mills    Level: developer
299a873a8cdSSam Reynolds 
300a873a8cdSSam Reynolds    Notes:
301a873a8cdSSam Reynolds     Each process has ALL the column reductions after the call. Because of the way this is computed each process gets all the values,
302a873a8cdSSam Reynolds     if each process wants only some of the values it should extract the ones it wants from the array.
303a873a8cdSSam Reynolds 
304a873a8cdSSam Reynolds   Developer Note:
305857cbf51SRichard Tran Mills     This routine is primarily intended as a back-end.
306857cbf51SRichard Tran Mills     MatGetColumnNorms(), MatGetColumnSums(), and MatGetColumnMeans() are implemented using this routine.
307a873a8cdSSam Reynolds 
308857cbf51SRichard Tran Mills .seealso: ReductionType, NormType, MatGetColumnNorms(), MatGetColumnSums(), MatGetColumnMeans()
309a873a8cdSSam Reynolds 
310a873a8cdSSam Reynolds @*/
311857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnReductions(Mat A,PetscInt type,PetscReal reductions[])
312a873a8cdSSam Reynolds {
313a873a8cdSSam Reynolds   PetscFunctionBegin;
314a873a8cdSSam Reynolds   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
315a873a8cdSSam Reynolds   if (A->ops->getcolumnreductions) {
3169566063dSJacob Faibussowitsch     PetscCall((*A->ops->getcolumnreductions)(A,type,reductions));
317a873a8cdSSam Reynolds   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not coded for this matrix type");
318a873a8cdSSam Reynolds   PetscFunctionReturn(0);
319a873a8cdSSam Reynolds }
320