Set of functions to manipulate saddle-point systems defined by blocks relying on the cs_cdo_system_t structures: matrix/vector multiplications. More...
#include "cs_defs.h"
#include <assert.h>
#include "bft_error.h"
#include "bft_mem.h"
#include "cs_array.h"
#include "cs_math.h"
#include "cs_saddle_system.h"
Functions | |
cs_real_t * | cs_saddle_system_b11_inv_diag (cs_lnum_t b11_max_size, cs_cdo_system_helper_t *sh) |
Retrieve the inverse of the diagonal of the (1,1)-block matrix The storage of a matrix is in a gather view and the resulting array is in a scatter view. More... | |
void | cs_saddle_system_b11_matvec (const cs_cdo_system_helper_t *sh, cs_real_t *vec, cs_real_t *matvec) |
Perform a matrix-vector multiplication for the (1,1) block when the input vector is in a scatter state. More... | |
void | cs_saddle_system_b12_matvec (const cs_cdo_system_helper_t *sh, const cs_real_t *x2, cs_real_t *res, bool reset_res) |
Compute the resulting vector of the operation m12*x2 The block (1,2) is stored in an unassembled way and corresponds to the (2,1) block. Thus, one considers a transposition of this block. x2 corresponds to a "scatter" view. More... | |
void | cs_saddle_system_b21_matvec (const cs_cdo_system_helper_t *sh, const cs_real_t *x1, cs_real_t *res) |
Compute the resulting vector of the operation m21*x1 The block (2,1) is stored in an unassembled way. x1 corresponds to a "scatter" view. More... | |
void | cs_saddle_system_matvec (const cs_cdo_system_helper_t *sh, cs_real_t *x1, cs_real_t *x2, cs_real_t *r1, cs_real_t *r2) |
Compute the matrix-vector operation for a saddle-point system r1 = M11.x1 + M12.x2 (result for the first set) r2 = M21.x1 (result for the second set) More... | |
void | cs_saddle_system_residual (const cs_cdo_system_helper_t *sh, cs_real_t *x1, cs_real_t *x2, cs_real_t *res1, cs_real_t *res2) |
Compute the residual of the saddle-point system res1 = rhs1 - M11.x1 - M12.x2 (residual for the first set) res2 = rhs2 - M21.x1 (residual for the second set) More... | |
Set of functions to manipulate saddle-point systems defined by blocks relying on the cs_cdo_system_t structures: matrix/vector multiplications.
cs_real_t* cs_saddle_system_b11_inv_diag | ( | cs_lnum_t | b11_max_size, |
cs_cdo_system_helper_t * | sh | ||
) |
Retrieve the inverse of the diagonal of the (1,1)-block matrix The storage of a matrix is in a gather view and the resulting array is in a scatter view.
[in] | b11_max_size | max size related to the (1,1) block |
[in] | sh | pointer to a system helper structure |
void cs_saddle_system_b11_matvec | ( | const cs_cdo_system_helper_t * | sh, |
cs_real_t * | vec, | ||
cs_real_t * | matvec | ||
) |
Perform a matrix-vector multiplication for the (1,1) block when the input vector is in a scatter state.
Thus, one performs a scatter --> gather (before the multiplication) and a gather --> scatter operation after the multiplication. One assumes that m11x1 is allocated to the right size. No check is done.
[in] | sh | pointer to a system helper structure |
[in,out] | vec | vector |
[in,out] | matvec | resulting vector for the matrix-vector product |
void cs_saddle_system_b12_matvec | ( | const cs_cdo_system_helper_t * | sh, |
const cs_real_t * | x2, | ||
cs_real_t * | res, | ||
bool | reset_res | ||
) |
Compute the resulting vector of the operation m12*x2 The block (1,2) is stored in an unassembled way and corresponds to the (2,1) block. Thus, one considers a transposition of this block. x2 corresponds to a "scatter" view.
[in] | sh | pointer to a system helper structure |
[in] | x2 | array to be multiplied |
[in,out] | res | result array storing the matrix.vector operation |
[in] | reset_res | false --> this is an update |
void cs_saddle_system_b21_matvec | ( | const cs_cdo_system_helper_t * | sh, |
const cs_real_t * | x1, | ||
cs_real_t * | res | ||
) |
Compute the resulting vector of the operation m21*x1 The block (2,1) is stored in an unassembled way. x1 corresponds to a "scatter" view.
[in] | sh | pointer to a system helper structure |
[in] | x1 | array to be multiplied |
[in,out] | res | result array storing the matrix.vector operation |
void cs_saddle_system_matvec | ( | const cs_cdo_system_helper_t * | sh, |
cs_real_t * | x1, | ||
cs_real_t * | x2, | ||
cs_real_t * | r1, | ||
cs_real_t * | r2 | ||
) |
Compute the matrix-vector operation for a saddle-point system r1 = M11.x1 + M12.x2 (result for the first set) r2 = M21.x1 (result for the second set)
[in] | sh | pointer to a system helper structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
[in,out] | r1 | result array for the first set of DoFs |
[in,out] | r2 | result array for the second set of DoFs |
void cs_saddle_system_residual | ( | const cs_cdo_system_helper_t * | sh, |
cs_real_t * | x1, | ||
cs_real_t * | x2, | ||
cs_real_t * | res1, | ||
cs_real_t * | res2 | ||
) |
Compute the residual of the saddle-point system res1 = rhs1 - M11.x1 - M12.x2 (residual for the first set) res2 = rhs2 - M21.x1 (residual for the second set)
[in] | sh | pointer to a system helper structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
[in,out] | res1 | resulting residual vector for the first set of DoFs |
[in,out] | res2 | resulting residual vector for the second set of DoFs |