xref: /petsc/src/mat/impls/centering/centering.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
1c1fff1c6SRichard Tran Mills 
2c1fff1c6SRichard Tran Mills #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
3c1fff1c6SRichard Tran Mills 
4d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Centering(Mat A, Vec xx, Vec yy)
5d71ae5a4SJacob Faibussowitsch {
6c1fff1c6SRichard Tran Mills   PetscScalar       *y;
7c1fff1c6SRichard Tran Mills   const PetscScalar *x;
8c1fff1c6SRichard Tran Mills   PetscScalar        sum, mean;
9c1fff1c6SRichard Tran Mills   PetscInt           i, m = A->rmap->n, size;
10c1fff1c6SRichard Tran Mills 
11c1fff1c6SRichard Tran Mills   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(VecSum(xx, &sum));
139566063dSJacob Faibussowitsch   PetscCall(VecGetSize(xx, &size));
14c1fff1c6SRichard Tran Mills   mean = sum / (PetscScalar)size;
159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
169566063dSJacob Faibussowitsch   PetscCall(VecGetArray(yy, &y));
17ad540459SPierre Jolivet   for (i = 0; i < m; i++) y[i] = x[i] - mean;
189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(yy, &y));
203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21c1fff1c6SRichard Tran Mills }
22c1fff1c6SRichard Tran Mills 
23c1fff1c6SRichard Tran Mills /*@
24c1fff1c6SRichard Tran Mills    MatCreateCentering - Creates a new matrix object that implements the (symmetric and idempotent) centering matrix, I - (1/N) * ones*ones'
25c1fff1c6SRichard Tran Mills 
2611a5261eSBarry Smith    Collective
27c1fff1c6SRichard Tran Mills 
28c1fff1c6SRichard Tran Mills    Input Parameters:
29c1fff1c6SRichard Tran Mills +  comm - MPI communicator
302ef1f0ffSBarry Smith .  n - number of local rows (or `PETSC_DECIDE` to have calculated if `N` is given)
31c1fff1c6SRichard Tran Mills            This value should be the same as the local size used in creating the
32c1fff1c6SRichard Tran Mills            y vector for the matrix-vector product y = Ax.
332ef1f0ffSBarry Smith -  N - number of global rows (or `PETSC_DETERMINE` to have calculated if `n` is given)
34c1fff1c6SRichard Tran Mills 
35c1fff1c6SRichard Tran Mills    Output Parameter:
36c1fff1c6SRichard Tran Mills .  C - the matrix
37c1fff1c6SRichard Tran Mills 
38c1fff1c6SRichard Tran Mills    Notes:
392ef1f0ffSBarry Smith    The entries of the matrix `C` are not explicitly stored. Instead, the new matrix
40c1fff1c6SRichard Tran Mills    object is a shell that simply performs the action of the centering matrix, i.e.,
41c1fff1c6SRichard Tran Mills    multiplying C*x subtracts the mean of the vector x from each of its elements.
42c1fff1c6SRichard Tran Mills    This is useful for preserving sparsity when mean-centering the columns of a
43c1fff1c6SRichard Tran Mills    matrix is required. For instance, to perform principal components analysis with
44c1fff1c6SRichard Tran Mills    a matrix A, the composite matrix C*A can be passed to a partial SVD solver.
45c1fff1c6SRichard Tran Mills 
46c1fff1c6SRichard Tran Mills    Level: intermediate
47c1fff1c6SRichard Tran Mills 
48*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateLRC()`, `MatCreateComposite()`
49c1fff1c6SRichard Tran Mills @*/
50d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateCentering(MPI_Comm comm, PetscInt n, PetscInt N, Mat *C)
51d71ae5a4SJacob Faibussowitsch {
52c1fff1c6SRichard Tran Mills   PetscMPIInt size;
53c1fff1c6SRichard Tran Mills 
54c1fff1c6SRichard Tran Mills   PetscFunctionBegin;
559566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, C));
569566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*C, n, n, N, N));
579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
589566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*C, MATCENTERING));
59c1fff1c6SRichard Tran Mills 
60c1fff1c6SRichard Tran Mills   (*C)->ops->mult                   = MatMult_Centering;
61c1fff1c6SRichard Tran Mills   (*C)->assembled                   = PETSC_TRUE;
62c1fff1c6SRichard Tran Mills   (*C)->preallocated                = PETSC_TRUE;
63b94d7dedSBarry Smith   (*C)->symmetric                   = PETSC_BOOL3_TRUE;
64b94d7dedSBarry Smith   (*C)->symmetry_eternal            = PETSC_TRUE;
65b94d7dedSBarry Smith   (*C)->structural_symmetry_eternal = PETSC_TRUE;
669566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*C));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68c1fff1c6SRichard Tran Mills }
69