xref: /petsc/src/mat/utils/getcolv.c (revision 857cbf51cdeafdb93f9c50bd4cde04c6075da3dc)
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;
31dfbe8321SBarry Smith   PetscErrorCode    ierr;
3238baddfdSBarry Smith   PetscInt          i,j,nz,N,Rs,Re,rs,re;
3338baddfdSBarry Smith   const PetscInt    *idx;
340925cdddSBarry Smith 
350925cdddSBarry Smith   PetscFunctionBegin;
360700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
370700a824SBarry Smith   PetscValidHeaderSpecific(yy,VEC_CLASSID,2);
389a971ab9SStefano Zampini   PetscValidLogicalCollectiveInt(A,col,3);
399a971ab9SStefano Zampini   if (col < 0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Requested negative column: %D",col);
400298fd71SBarry Smith   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
419a971ab9SStefano Zampini   if (col >= N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Requested column %D larger than number columns in matrix %D",col,N);
4282bf6240SBarry Smith   ierr = MatGetOwnershipRange(A,&Rs,&Re);CHKERRQ(ierr);
4382bf6240SBarry Smith   ierr = VecGetOwnershipRange(yy,&rs,&re);CHKERRQ(ierr);
44e32f2f54SBarry Smith   if (Rs != rs || Re != re) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Matrix %D %D does not have same ownership range (size) as vector %D %D",Rs,Re,rs,re);
4582bf6240SBarry Smith 
468d0534beSBarry Smith   if (A->ops->getcolumnvector) {
478d0534beSBarry Smith     ierr = (*A->ops->getcolumnvector)(A,yy,col);CHKERRQ(ierr);
488d0534beSBarry Smith   } else {
498aa348c1SBarry Smith     ierr = VecSet(yy,0.0);CHKERRQ(ierr);
5082bf6240SBarry Smith     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
519a971ab9SStefano Zampini     /* TODO for general matrices */
5282bf6240SBarry Smith     for (i=Rs; i<Re; i++) {
5382bf6240SBarry Smith       ierr = MatGetRow(A,i,&nz,&idx,&v);CHKERRQ(ierr);
5482bf6240SBarry Smith       if (nz && idx[0] <= col) {
5582bf6240SBarry Smith         /*
5682bf6240SBarry Smith           Should use faster search here
5782bf6240SBarry Smith         */
5882bf6240SBarry Smith         for (j=0; j<nz; j++) {
5982bf6240SBarry Smith           if (idx[j] >= col) {
6082bf6240SBarry Smith             if (idx[j] == col) y[i-rs] = v[j];
6182bf6240SBarry Smith             break;
620925cdddSBarry Smith           }
630925cdddSBarry Smith         }
640925cdddSBarry Smith       }
6582bf6240SBarry Smith       ierr = MatRestoreRow(A,i,&nz,&idx,&v);CHKERRQ(ierr);
660925cdddSBarry Smith     }
6782bf6240SBarry Smith     ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
688d0534beSBarry Smith   }
690925cdddSBarry Smith   PetscFunctionReturn(0);
700925cdddSBarry Smith }
71242f1d38SBarry Smith 
7286819fdcSBarry Smith /*@
730716a85fSBarry Smith    MatGetColumnNorms - Gets the norms of each column of a sparse or dense matrix.
74242f1d38SBarry Smith 
75242f1d38SBarry Smith    Input Parameter:
76242f1d38SBarry Smith +  A - the matrix
77242f1d38SBarry Smith -  type - NORM_2, NORM_1 or NORM_INFINITY
78242f1d38SBarry Smith 
79242f1d38SBarry Smith    Output Parameter:
80242f1d38SBarry Smith .  norms - an array as large as the TOTAL number of columns in the matrix
81242f1d38SBarry Smith 
82f6680f47SSatish Balay    Level: intermediate
83f6680f47SSatish Balay 
8495452b02SPatrick Sanan    Notes:
8595452b02SPatrick Sanan     Each process has ALL the column norms after the call. Because of the way this is computed each process gets all the values,
8686819fdcSBarry Smith     if each process wants only some of the values it should extract the ones it wants from the array.
8786819fdcSBarry Smith 
8840b1048aSBarry Smith .seealso: NormType, MatNorm()
8986819fdcSBarry Smith 
9086819fdcSBarry Smith @*/
916b9dab40SJed Brown PetscErrorCode MatGetColumnNorms(Mat A,NormType type,PetscReal norms[])
92242f1d38SBarry Smith {
93*857cbf51SRichard Tran Mills   /* NOTE: MatGetColumnNorms() could simply be a macro that calls MatGetColumnReductions().
94*857cbf51SRichard Tran Mills    * I've kept this as a function because it allows slightly more in the way of error checking,
95*857cbf51SRichard Tran Mills    * erroring out if MatGetColumnNorms() is not called with a valid NormType. */
96242f1d38SBarry Smith   PetscErrorCode ierr;
97242f1d38SBarry Smith 
98242f1d38SBarry Smith   PetscFunctionBegin;
99*857cbf51SRichard Tran Mills   if (type == NORM_2 || type == NORM_1 || type == NORM_FROBENIUS || type == NORM_INFINITY || type == NORM_1_AND_2) {
100*857cbf51SRichard Tran Mills     ierr = MatGetColumnReductions(A,type,norms);CHKERRQ(ierr);
101*857cbf51SRichard Tran Mills   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType");
102*857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
103a873a8cdSSam Reynolds }
104*857cbf51SRichard Tran Mills 
105*857cbf51SRichard Tran Mills /*@
106*857cbf51SRichard Tran Mills    MatGetColumnSumsRealPart - Gets the sums of the real part of each column of a sparse or dense matrix.
107*857cbf51SRichard Tran Mills 
108*857cbf51SRichard Tran Mills    Input Parameter:
109*857cbf51SRichard Tran Mills .  A - the matrix
110*857cbf51SRichard Tran Mills 
111*857cbf51SRichard Tran Mills    Output Parameter:
112*857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
113*857cbf51SRichard Tran Mills 
114*857cbf51SRichard Tran Mills    Level: intermediate
115*857cbf51SRichard Tran Mills 
116*857cbf51SRichard Tran Mills    Notes:
117*857cbf51SRichard 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,
118*857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
119*857cbf51SRichard Tran Mills 
120*857cbf51SRichard Tran Mills .seealso: MatGetColumnSumsImaginaryPart(), VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions()
121*857cbf51SRichard Tran Mills 
122*857cbf51SRichard Tran Mills @*/
123*857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSumsRealPart(Mat A,PetscReal sums[])
124*857cbf51SRichard Tran Mills {
125*857cbf51SRichard Tran Mills   PetscErrorCode ierr;
126*857cbf51SRichard Tran Mills 
127*857cbf51SRichard Tran Mills   PetscFunctionBegin;
128*857cbf51SRichard Tran Mills   ierr = MatGetColumnReductions(A,REDUCTION_SUM_REALPART,sums);CHKERRQ(ierr);
129*857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
130*857cbf51SRichard Tran Mills }
131*857cbf51SRichard Tran Mills 
132*857cbf51SRichard Tran Mills /*@
133*857cbf51SRichard Tran Mills    MatGetColumnSumsImaginaryPart - Gets the sums of the imaginary part of each column of a sparse or dense matrix.
134*857cbf51SRichard Tran Mills 
135*857cbf51SRichard Tran Mills    Input Parameter:
136*857cbf51SRichard Tran Mills .  A - the matrix
137*857cbf51SRichard Tran Mills 
138*857cbf51SRichard Tran Mills    Output Parameter:
139*857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
140*857cbf51SRichard Tran Mills 
141*857cbf51SRichard Tran Mills    Level: intermediate
142*857cbf51SRichard Tran Mills 
143*857cbf51SRichard Tran Mills    Notes:
144*857cbf51SRichard 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,
145*857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
146*857cbf51SRichard Tran Mills 
147*857cbf51SRichard Tran Mills .seealso: MatGetColumnSumsRealPart(), VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions()
148*857cbf51SRichard Tran Mills 
149*857cbf51SRichard Tran Mills @*/
150*857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSumsImaginaryPart(Mat A,PetscReal sums[])
151*857cbf51SRichard Tran Mills {
152*857cbf51SRichard Tran Mills   PetscErrorCode ierr;
153*857cbf51SRichard Tran Mills 
154*857cbf51SRichard Tran Mills   PetscFunctionBegin;
155*857cbf51SRichard Tran Mills   ierr = MatGetColumnReductions(A,REDUCTION_SUM_IMAGINARYPART,sums);CHKERRQ(ierr);
156*857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
157*857cbf51SRichard Tran Mills }
158*857cbf51SRichard Tran Mills 
159*857cbf51SRichard Tran Mills /*@
160*857cbf51SRichard Tran Mills    MatGetColumnSums - Gets the sums of each column of a sparse or dense matrix.
161*857cbf51SRichard Tran Mills 
162*857cbf51SRichard Tran Mills    Input Parameter:
163*857cbf51SRichard Tran Mills .  A - the matrix
164*857cbf51SRichard Tran Mills 
165*857cbf51SRichard Tran Mills    Output Parameter:
166*857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
167*857cbf51SRichard Tran Mills 
168*857cbf51SRichard Tran Mills    Level: intermediate
169*857cbf51SRichard Tran Mills 
170*857cbf51SRichard Tran Mills    Notes:
171*857cbf51SRichard 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,
172*857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
173*857cbf51SRichard Tran Mills 
174*857cbf51SRichard Tran Mills .seealso: VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions()
175*857cbf51SRichard Tran Mills 
176*857cbf51SRichard Tran Mills @*/
177*857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSums(Mat A,PetscScalar sums[])
178*857cbf51SRichard Tran Mills {
179*857cbf51SRichard Tran Mills   PetscErrorCode ierr;
180*857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
181*857cbf51SRichard Tran Mills   PetscInt       i,n;
182*857cbf51SRichard Tran Mills   PetscReal      *work;
183*857cbf51SRichard Tran Mills #endif
184*857cbf51SRichard Tran Mills 
185*857cbf51SRichard Tran Mills   PetscFunctionBegin;
186*857cbf51SRichard Tran Mills 
187*857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
188*857cbf51SRichard Tran Mills   ierr = MatGetColumnSumsRealPart(A,sums);CHKERRQ(ierr);
189*857cbf51SRichard Tran Mills #else
190*857cbf51SRichard Tran Mills   ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr);
191*857cbf51SRichard Tran Mills   ierr = PetscArrayzero(sums,n);CHKERRQ(ierr);
192*857cbf51SRichard Tran Mills   ierr = PetscCalloc1(n,&work);CHKERRQ(ierr);
193*857cbf51SRichard Tran Mills   ierr = MatGetColumnSumsRealPart(A,work);CHKERRQ(ierr);
194*857cbf51SRichard Tran Mills   for (i=0; i<n; i++) sums[i] = work[i];
195*857cbf51SRichard Tran Mills   ierr = MatGetColumnSumsImaginaryPart(A,work);CHKERRQ(ierr);
196*857cbf51SRichard Tran Mills   for (i=0; i<n; i++) sums[i] += work[i]*PETSC_i;
197*857cbf51SRichard Tran Mills   ierr = PetscFree(work);CHKERRQ(ierr);
198*857cbf51SRichard Tran Mills #endif
199*857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
200*857cbf51SRichard Tran Mills }
201*857cbf51SRichard Tran Mills 
202*857cbf51SRichard Tran Mills /*@
203*857cbf51SRichard Tran Mills    MatGetColumnMeansRealPart - Gets the arithmetic means of the real part of each column of a sparse or dense matrix.
204*857cbf51SRichard Tran Mills 
205*857cbf51SRichard Tran Mills    Input Parameter:
206*857cbf51SRichard Tran Mills .  A - the matrix
207*857cbf51SRichard Tran Mills 
208*857cbf51SRichard Tran Mills    Output Parameter:
209*857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
210*857cbf51SRichard Tran Mills 
211*857cbf51SRichard Tran Mills    Level: intermediate
212*857cbf51SRichard Tran Mills 
213*857cbf51SRichard Tran Mills    Notes:
214*857cbf51SRichard 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,
215*857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
216*857cbf51SRichard Tran Mills 
217*857cbf51SRichard Tran Mills .seealso: MatGetColumnMeansImaginaryPart(), VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions()
218*857cbf51SRichard Tran Mills 
219*857cbf51SRichard Tran Mills @*/
220*857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeansRealPart(Mat A,PetscReal means[])
221*857cbf51SRichard Tran Mills {
222*857cbf51SRichard Tran Mills   PetscErrorCode ierr;
223*857cbf51SRichard Tran Mills 
224*857cbf51SRichard Tran Mills   PetscFunctionBegin;
225*857cbf51SRichard Tran Mills   ierr = MatGetColumnReductions(A,REDUCTION_MEAN_REALPART,means);CHKERRQ(ierr);
226*857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
227*857cbf51SRichard Tran Mills }
228*857cbf51SRichard Tran Mills 
229*857cbf51SRichard Tran Mills /*@
230*857cbf51SRichard Tran Mills    MatGetColumnMeansImaginaryPart - Gets the arithmetic means of the imaginary part of each column of a sparse or dense matrix.
231*857cbf51SRichard Tran Mills 
232*857cbf51SRichard Tran Mills    Input Parameter:
233*857cbf51SRichard Tran Mills .  A - the matrix
234*857cbf51SRichard Tran Mills 
235*857cbf51SRichard Tran Mills    Output Parameter:
236*857cbf51SRichard Tran Mills .  sums - an array as large as the TOTAL number of columns in the matrix
237*857cbf51SRichard Tran Mills 
238*857cbf51SRichard Tran Mills    Level: intermediate
239*857cbf51SRichard Tran Mills 
240*857cbf51SRichard Tran Mills    Notes:
241*857cbf51SRichard 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,
242*857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
243*857cbf51SRichard Tran Mills 
244*857cbf51SRichard Tran Mills .seealso: MatGetColumnMeansRealPart(), VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions()
245*857cbf51SRichard Tran Mills 
246*857cbf51SRichard Tran Mills @*/
247*857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeansImaginaryPart(Mat A,PetscReal means[])
248*857cbf51SRichard Tran Mills {
249*857cbf51SRichard Tran Mills   PetscErrorCode ierr;
250*857cbf51SRichard Tran Mills 
251*857cbf51SRichard Tran Mills   PetscFunctionBegin;
252*857cbf51SRichard Tran Mills   ierr = MatGetColumnReductions(A,REDUCTION_MEAN_IMAGINARYPART,means);CHKERRQ(ierr);
253*857cbf51SRichard Tran Mills   PetscFunctionReturn(0);
254*857cbf51SRichard Tran Mills }
255*857cbf51SRichard Tran Mills 
256*857cbf51SRichard Tran Mills /*@
257*857cbf51SRichard Tran Mills    MatGetColumnMeans - Gets the arithmetic means of each column of a sparse or dense matrix.
258*857cbf51SRichard Tran Mills 
259*857cbf51SRichard Tran Mills    Input Parameter:
260*857cbf51SRichard Tran Mills .  A - the matrix
261*857cbf51SRichard Tran Mills 
262*857cbf51SRichard Tran Mills    Output Parameter:
263*857cbf51SRichard Tran Mills .  means - an array as large as the TOTAL number of columns in the matrix
264*857cbf51SRichard Tran Mills 
265*857cbf51SRichard Tran Mills    Level: intermediate
266*857cbf51SRichard Tran Mills 
267*857cbf51SRichard Tran Mills    Notes:
268*857cbf51SRichard 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,
269*857cbf51SRichard Tran Mills     if each process wants only some of the values it should extract the ones it wants from the array.
270*857cbf51SRichard Tran Mills 
271*857cbf51SRichard Tran Mills .seealso: VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions()
272*857cbf51SRichard Tran Mills 
273*857cbf51SRichard Tran Mills @*/
274*857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeans(Mat A,PetscScalar means[])
275*857cbf51SRichard Tran Mills {
276*857cbf51SRichard Tran Mills   PetscErrorCode ierr;
277*857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX)
278*857cbf51SRichard Tran Mills   PetscInt       i,n;
279*857cbf51SRichard Tran Mills   PetscReal      *work;
280*857cbf51SRichard Tran Mills #endif
281*857cbf51SRichard Tran Mills 
282*857cbf51SRichard Tran Mills   PetscFunctionBegin;
283*857cbf51SRichard Tran Mills 
284*857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX)
285*857cbf51SRichard Tran Mills   ierr = MatGetColumnMeansRealPart(A,means);CHKERRQ(ierr);
286*857cbf51SRichard Tran Mills #else
287*857cbf51SRichard Tran Mills   ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr);
288*857cbf51SRichard Tran Mills   ierr = PetscArrayzero(means,n);CHKERRQ(ierr);
289*857cbf51SRichard Tran Mills   ierr = PetscCalloc1(n,&work);CHKERRQ(ierr);
290*857cbf51SRichard Tran Mills   ierr = MatGetColumnMeansRealPart(A,work);CHKERRQ(ierr);
291*857cbf51SRichard Tran Mills   for (i=0; i<n; i++) means[i] = work[i];
292*857cbf51SRichard Tran Mills   ierr = MatGetColumnMeansImaginaryPart(A,work);CHKERRQ(ierr);
293*857cbf51SRichard Tran Mills   for (i=0; i<n; i++) means[i] += work[i]*PETSC_i;
294*857cbf51SRichard Tran Mills   ierr = PetscFree(work);CHKERRQ(ierr);
295*857cbf51SRichard Tran Mills #endif
296242f1d38SBarry Smith   PetscFunctionReturn(0);
297242f1d38SBarry Smith }
298242f1d38SBarry Smith 
299a873a8cdSSam Reynolds /*@
300a873a8cdSSam Reynolds     MatGetColumnReductions - Gets the reductions of each column of a sparse or dense matrix.
301a873a8cdSSam Reynolds 
302a873a8cdSSam Reynolds   Input Parameter:
303a873a8cdSSam Reynolds +  A - the matrix
304*857cbf51SRichard Tran Mills -  type - A constant defined in NormType or ReductionType: NORM_2, NORM_1, NORM_INFINITY, REDUCTION_SUM_REALPART,
305*857cbf51SRichard Tran Mills           REDUCTION_SUM_IMAGINARYPART, REDUCTION_MEAN_REALPART, REDUCTION_MEAN_IMAGINARYPART
306a873a8cdSSam Reynolds 
307a873a8cdSSam Reynolds   Output Parameter:
308a873a8cdSSam Reynolds .  reductions - an array as large as the TOTAL number of columns in the matrix
309a873a8cdSSam Reynolds 
310*857cbf51SRichard Tran Mills    Level: developer
311a873a8cdSSam Reynolds 
312a873a8cdSSam Reynolds    Notes:
313a873a8cdSSam Reynolds     Each process has ALL the column reductions after the call. Because of the way this is computed each process gets all the values,
314a873a8cdSSam Reynolds     if each process wants only some of the values it should extract the ones it wants from the array.
315a873a8cdSSam Reynolds 
316a873a8cdSSam Reynolds   Developer Note:
317*857cbf51SRichard Tran Mills     This routine is primarily intended as a back-end.
318*857cbf51SRichard Tran Mills     MatGetColumnNorms(), MatGetColumnSums(), and MatGetColumnMeans() are implemented using this routine.
319a873a8cdSSam Reynolds 
320*857cbf51SRichard Tran Mills .seealso: ReductionType, NormType, MatGetColumnNorms(), MatGetColumnSums(), MatGetColumnMeans()
321a873a8cdSSam Reynolds 
322a873a8cdSSam Reynolds @*/
323*857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnReductions(Mat A,PetscInt type,PetscReal reductions[])
324a873a8cdSSam Reynolds {
325a873a8cdSSam Reynolds   PetscErrorCode ierr;
326a873a8cdSSam Reynolds 
327a873a8cdSSam Reynolds   PetscFunctionBegin;
328a873a8cdSSam Reynolds   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
329a873a8cdSSam Reynolds   if (A->ops->getcolumnreductions) {
330a873a8cdSSam Reynolds     ierr = (*A->ops->getcolumnreductions)(A,type,reductions);CHKERRQ(ierr);
331a873a8cdSSam Reynolds   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not coded for this matrix type");
332a873a8cdSSam Reynolds   PetscFunctionReturn(0);
333a873a8cdSSam Reynolds }
334