9.0
general documentation
Loading...
Searching...
No Matches
cs_saddle_system.h File Reference
#include "base/cs_defs.h"
#include "cdo/cs_cdo_system.h"
Include dependency graph for cs_saddle_system.h:

Go to the source code of this file.

Functions

cs_real_tcs_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.
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.
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.
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.
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)
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)

Function Documentation

◆ cs_saddle_system_b11_inv_diag()

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.

Parameters
[in]b11_max_sizemax size related to the (1,1) block
[in]shpointer to a system helper structure
Returns
a pointer to the computed array (scatter view)

◆ cs_saddle_system_b11_matvec()

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.

Parameters
[in]shpointer to a system helper structure
[in,out]vecvector
[in,out]matvecresulting vector for the matrix-vector product

◆ cs_saddle_system_b12_matvec()

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.

Parameters
[in]shpointer to a system helper structure
[in]x2array to be multiplied
[in,out]resresult array storing the matrix.vector operation
[in]reset_resfalse --> this is an update

◆ cs_saddle_system_b21_matvec()

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.

Parameters
[in]shpointer to a system helper structure
[in]x1array to be multiplied
[in,out]resresult array storing the matrix.vector operation

◆ cs_saddle_system_matvec()

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)

Parameters
[in]shpointer to a system helper structure
[in,out]x1array for the first part
[in,out]x2array for the second part
[in,out]r1result array for the first set of DoFs
[in,out]r2result array for the second set of DoFs

◆ cs_saddle_system_residual()

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)

Parameters
[in]shpointer to a system helper structure
[in,out]x1array for the first part
[in,out]x2array for the second part
[in,out]res1resulting residual vector for the first set of DoFs
[in,out]res2resulting residual vector for the second set of DoFs