:orphan: # MatGetColumnReductions Gets the reductions of each column of a sparse or dense matrix. ## Synopsis ``` #include "petscmat.h" PetscErrorCode MatGetColumnReductions(Mat A, PetscInt type, PetscReal reductions[]) ``` ## Input Parameters - ***A -*** the matrix - ***type -*** A constant defined in `NormType` or `ReductionType`: `NORM_2`, `NORM_1`, `NORM_INFINITY`, `REDUCTION_SUM_REALPART`, `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_MEAN_REALPART`, `REDUCTION_MEAN_IMAGINARYPART` ## Output Parameter - ***reductions -*** an array as large as the TOTAL number of columns in the matrix ## Note Each process has ALL the column reductions after the call. Because of the way this is computed each process gets all the values, if each process wants only some of the values it should extract the ones it wants from the array. ## Developer Note This routine is primarily intended as a back-end. `MatGetColumnNorms()`, `MatGetColumnSums()`, and `MatGetColumnMeans()` are implemented using this routine. ## See Also [](ch_matrices), `Mat`, `ReductionType`, `NormType`, `MatGetColumnNorms()`, `MatGetColumnSums()`, `MatGetColumnMeans()` ## Level developer ## Location src/mat/utils/getcolv.c ## Implementations MatGetColumnReductions_MPIAIJ in src/mat/impls/aij/mpi/mpiaij.c
MatGetColumnReductions_SeqAIJ in src/mat/impls/aij/seq/aij.c
MatGetColumnReductions_MPIBAIJ in src/mat/impls/baij/mpi/mpibaij.c
MatGetColumnReductions_SeqBAIJ in src/mat/impls/baij/seq/baij.c
MatGetColumnReductions_MPIDense in src/mat/impls/dense/mpi/mpidense.c
MatGetColumnReductions_SeqDense in src/mat/impls/dense/seq/dense.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/utils/getcolv.c) [Index of all Mat routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)