#include "cs_defs.h"
#include "cs_cdo_blas.h"
#include "cs_cdo_system.h"
#include "cs_iter_algo.h"
#include "cs_param_saddle.h"
#include "cs_sles.h"
Go to the source code of this file.
Data Structures | |
struct | cs_saddle_solver_t |
struct | cs_saddle_solver_context_alu_t |
struct | cs_saddle_solver_context_block_pcd_t |
struct | cs_saddle_solver_context_gkb_t |
struct | cs_saddle_solver_context_notay_t |
struct | cs_saddle_solver_context_uzawa_cg_t |
struct | cs_saddle_solver_context_simple_t |
Typedefs | |
typedef void() | cs_saddle_solver_matvec_t(cs_lnum_t n2_dofs, const cs_real_t *vec, const cs_adjacency_t *mat_adj, const cs_real_t *mat_op, cs_real_t *matvec) |
Generic function prototype to perform a matrix vector operation This operation takes place between an unassembled matrix and a vector. More... | |
Functions | |
int | cs_saddle_solver_get_n_systems (void) |
Retrieve the number of saddle-point systems which have been added. More... | |
cs_saddle_solver_t * | cs_saddle_solver_by_id (int id) |
Get a pointer to saddle-point solver from its id. More... | |
cs_saddle_solver_t * | cs_saddle_solver_add (cs_lnum_t n1_elts, int n1_dofs_by_elt, cs_lnum_t n2_elts, int n2_dofs_by_elt, const cs_param_saddle_t *saddlep, cs_cdo_system_helper_t *sh, cs_sles_t *main_sles) |
Add a new solver for solving a saddle-point problem. More... | |
void | cs_saddle_solver_finalize (void) |
Free all remaining structures related to saddle-point solvers. More... | |
void | cs_saddle_solver_free (cs_saddle_solver_t **p_solver) |
Free a cs_saddle_solver_t structure. More... | |
void | cs_saddle_solver_clean (cs_saddle_solver_t *solver) |
Free/reset only a part of a cs_saddle_solver_t structure. More... | |
void | cs_saddle_solver_update_monitoring (cs_saddle_solver_t *solver, unsigned n_iter) |
Update the current monitoring state with n_iter. More... | |
void | cs_saddle_solver_log_monitoring (void) |
Log the monitoring performance of all saddle-point systems. More... | |
cs_real_t * | cs_saddle_solver_m11_inv_lumped (cs_saddle_solver_t *solver, const cs_matrix_t *m11, const cs_range_set_t *b11_rset, cs_sles_t *xtra_sles, int *n_iter) |
Retrieve the lumped matrix 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 scatter view. More... | |
void | cs_saddle_solver_m12_multiply_vector (cs_lnum_t n2_dofs, const cs_real_t *x2, const cs_adjacency_t *m21_adj, const cs_real_t *m21_val, cs_real_t *m12x2) |
Compute the resulting vector of the operation m12*x2 The stride is equal to 3 for the operator m21 (unassembled) x2 corresponds to a "scatter" view. More... | |
void | cs_saddle_solver_m12_multiply_scalar (cs_lnum_t n2_elts, const cs_real_t *x2, const cs_adjacency_t *m21_adj, const cs_real_t *m21_val, cs_real_t *m12x2) |
Compute the resulting vector of the operation m12*x2 The stride is equal to 1 for the operator m21 (unassembled) x2 corresponds to a "scatter" view. More... | |
void | cs_saddle_solver_m21_multiply_vector (cs_lnum_t n2_dofs, const cs_real_t *x1, const cs_adjacency_t *m21_adj, const cs_real_t *m21_val, cs_real_t *m21x1) |
Compute the resulting vector of the operation m21*x1 The stride is equal to 3 for the operator m21 operator x1 corresponds to a "scatter" view. More... | |
void | cs_saddle_solver_m21_multiply_scalar (cs_lnum_t n2_dofs, const cs_real_t *x1, const cs_adjacency_t *m21_adj, const cs_real_t *m21_val, cs_real_t *m21x1) |
Compute the resulting vector of the operation m21*x1 The stride is equal to 1 for the operator m21 operator x1 corresponds to a "scatter" view. More... | |
void | cs_saddle_solver_context_alu_create (cs_saddle_solver_t *solver) |
Create and initialize the context structure for an algorithm related to the ALU algorithm. More... | |
void | cs_saddle_solver_context_alu_clean (cs_saddle_solver_context_alu_t *ctx) |
Free main memory consuming part of the context structure associated to an ALU algorithm. More... | |
void | cs_saddle_solver_context_alu_free (cs_saddle_solver_context_alu_t **p_ctx) |
Free the context structure associated to an ALU algorithm. More... | |
void | cs_saddle_solver_context_block_pcd_create (cs_lnum_t b22_max_size, cs_saddle_solver_t *solver) |
Create and initialize the context structure for a block preconditioner used in combination with a Krylov solver such as MINRES or GCR. More... | |
void | cs_saddle_solver_context_block_pcd_clean (cs_saddle_solver_context_block_pcd_t *ctx) |
Free the context structure for a block preconditioner used in combination with a Krylov solver such as MINRES or GCR. More... | |
void | cs_saddle_solver_context_block_pcd_free (cs_saddle_solver_context_block_pcd_t **p_ctx) |
Free the context structure for a block preconditioner used in combination with a Krylov solver such as MINRES or GCR. More... | |
void | cs_saddle_solver_context_gkb_create (cs_saddle_solver_t *solver) |
Create and initialize the context structure for the GKB algorithm. More... | |
void | cs_saddle_solver_context_gkb_clean (cs_saddle_solver_context_gkb_t *ctx) |
Free the main memory consuming part of the context structure associated to the GKB algorithm. More... | |
void | cs_saddle_solver_context_gkb_free (cs_saddle_solver_context_gkb_t **p_ctx) |
Free the context structure associated to the GKB algorithm. More... | |
void | cs_saddle_solver_context_notay_create (cs_saddle_solver_t *solver) |
Create and initialize the context structure for the algorithm relying on the Notay's algebraic transformation. More... | |
void | cs_saddle_solver_context_uzawa_cg_create (cs_lnum_t b22_max_size, cs_saddle_solver_t *solver) |
Create and initialize the context structure for an algorithm related to the Uzawa-CG algorithm. More... | |
void | cs_saddle_solver_context_uzawa_cg_clean (cs_saddle_solver_context_uzawa_cg_t *ctx) |
Free main memory consuming part of the context structure associated to a Uzawa-CG algorithm. More... | |
void | cs_saddle_solver_context_uzawa_cg_free (cs_saddle_solver_context_uzawa_cg_t **p_ctx) |
Free the context structure associated to a Uzawa-CG algorithm. More... | |
void | cs_saddle_solver_context_simple_create (cs_lnum_t b22_max_size, cs_saddle_solver_t *solver) |
Create and initialize the context structure for an algorithm related to the SIMPLE algorithm. More... | |
void | cs_saddle_solver_context_simple_clean (cs_saddle_solver_context_simple_t *ctx) |
Free main memory consuming part of the context structure associated to a SIMPLE algorithm. More... | |
void | cs_saddle_solver_context_simple_free (cs_saddle_solver_context_simple_t **p_ctx) |
Free the context structure associated to a Uzawa-CG algorithm. More... | |
void | cs_saddle_solver_alu_incr (cs_saddle_solver_t *solver, cs_real_t *x1, cs_real_t *x2) |
Apply the Augmented Lagrangian-Uzawa algorithm to a saddle point problem (the system is stored in a hybrid way). The stride is equal to 1 for the matrix (db_size[3] = 1) and the vector. More... | |
void | cs_saddle_solver_gcr (cs_saddle_solver_t *solver, cs_real_t *x1, cs_real_t *x2) |
Apply the GCR algorithm to a saddle point problem (the system is stored in a hybrid way). The stride is equal to 1 for the matrix in the (1,1) block (db_size[3] = 1). More... | |
void | cs_saddle_solver_gkb_inhouse (cs_saddle_solver_t *solver, cs_real_t *x1, cs_real_t *x2) |
Apply the GKB algorithm to a saddle point problem (the system is stored in a hybrid way). The stride is equal to 1 for the matrix in the (1,1) block (db_size[3] = 1). More... | |
void | cs_saddle_solver_sles_full_system (cs_saddle_solver_t *solver, cs_real_t *x1, cs_real_t *x2) |
Apply an (external) solver to solve a saddle point problem (the system is stored in a monolithic way). More... | |
void | cs_saddle_solver_minres (cs_saddle_solver_t *solver, cs_real_t *x1, cs_real_t *x2) |
Apply the MINRES algorithm to a saddle point problem (the system is stored in a hybrid way). The stride is equal to 1 for the matrix (db_size[3] = 1) and the vector. More... | |
void | cs_saddle_solver_notay (cs_saddle_solver_t *solver, cs_real_t *x1, cs_real_t *x2) |
Apply the Notay's transformation algorithm to solve a saddle point problem (the system is stored in a monolithic way). More... | |
void | cs_saddle_solver_uzawa_cg (cs_saddle_solver_t *solver, cs_real_t *x1, cs_real_t *x2) |
Apply the Uzawa-CG algorithm to solve a saddle point problem (the system is stored in a hybrid way). This algorithm is based on Koko's paper "Uzawa conjugate gradient method for the Stokes problem: Matlab
implementation with P1-iso-P2/ P1 finite element". More... | |
void | cs_saddle_solver_simple (cs_saddle_solver_t *solver, cs_real_t *x1, cs_real_t *x2) |
Apply the SIMPLE-like algorithm to solve a saddle point problem. More... | |
typedef void() cs_saddle_solver_matvec_t(cs_lnum_t n2_dofs, const cs_real_t *vec, const cs_adjacency_t *mat_adj, const cs_real_t *mat_op, cs_real_t *matvec) |
Generic function prototype to perform a matrix vector operation This operation takes place between an unassembled matrix and a vector.
[in] | n2_dofs | number of (scatter) DoFs for the (2,2)-block |
[in] | vec | array of values |
[in] | mat_adj | adjacency related to the matrix operator |
[in] | mat_val | values associated to the matrix operator |
[in,out] | matvec | resulting vector (have to be allocated) |
cs_saddle_solver_t * cs_saddle_solver_add | ( | cs_lnum_t | n1_elts, |
int | n1_dofs_by_elt, | ||
cs_lnum_t | n2_elts, | ||
int | n2_dofs_by_elt, | ||
const cs_param_saddle_t * | saddlep, | ||
cs_cdo_system_helper_t * | sh, | ||
cs_sles_t * | main_sles | ||
) |
Add a new solver for solving a saddle-point problem.
[in] | n1_elts | number of elements associated to the (1,1)-block |
[in] | n1_dofs_by_elt | number of DoFs by elements in the (1,1)-block |
[in] | n2_elts | number of elements associated to the (2,2)-block |
[in] | n2_dofs_by_elt | number of DoFs by elements in the (2,2)-block |
[in] | saddlep | set of parameters for the saddle-point solver |
[in] | sh | pointer to a system helper structure |
[in] | main_sles | pointer to the main SLES structure related to this saddle-point problem |
void cs_saddle_solver_alu_incr | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | x1, | ||
cs_real_t * | x2 | ||
) |
Apply the Augmented Lagrangian-Uzawa algorithm to a saddle point problem (the system is stored in a hybrid way). The stride is equal to 1 for the matrix (db_size[3] = 1) and the vector.
[in,out] | solver | pointer to a saddle_point solver structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
cs_saddle_solver_t * cs_saddle_solver_by_id | ( | int | id | ) |
Get a pointer to saddle-point solver from its id.
[in] | id | id of the saddle-point system |
Get a pointer to saddle-point solver from its id.
[in] | id | id of the saddle-point system |
void cs_saddle_solver_clean | ( | cs_saddle_solver_t * | solver | ) |
Free/reset only a part of a cs_saddle_solver_t structure.
[in,out] | solver | double pointer to the structure to free |
void cs_saddle_solver_context_alu_clean | ( | cs_saddle_solver_context_alu_t * | ctx | ) |
Free main memory consuming part of the context structure associated to an ALU algorithm.
[in,out] | ctx | pointer to the context structure to clean |
void cs_saddle_solver_context_alu_create | ( | cs_saddle_solver_t * | solver | ) |
Create and initialize the context structure for an algorithm related to the ALU algorithm.
[in,out] | solver | pointer to a saddle-point solver structure |
void cs_saddle_solver_context_alu_free | ( | cs_saddle_solver_context_alu_t ** | p_ctx | ) |
Free the context structure associated to an ALU algorithm.
[in,out] | p_ctx | double pointer to the context structure to free |
void cs_saddle_solver_context_block_pcd_clean | ( | cs_saddle_solver_context_block_pcd_t * | ctx | ) |
Free the context structure for a block preconditioner used in combination with a Krylov solver such as MINRES or GCR.
[in,out] | ctx | pointer to the context structure to clean |
Free the context structure for a block preconditioner used in combination with a Krylov solver such as MINRES or GCR.
[in,out] | ctx | pointer to the context structure to free |
void cs_saddle_solver_context_block_pcd_create | ( | cs_lnum_t | b22_max_size, |
cs_saddle_solver_t * | solver | ||
) |
Create and initialize the context structure for a block preconditioner used in combination with a Krylov solver such as MINRES or GCR.
[in] | b22_max_size | max. size for the second part of unknows |
[in,out] | solver | pointer to a saddle-point solver structure |
void cs_saddle_solver_context_block_pcd_free | ( | cs_saddle_solver_context_block_pcd_t ** | p_ctx | ) |
Free the context structure for a block preconditioner used in combination with a Krylov solver such as MINRES or GCR.
[in,out] | p_ctx | double pointer to the context structure to free |
void cs_saddle_solver_context_gkb_clean | ( | cs_saddle_solver_context_gkb_t * | ctx | ) |
Free the main memory consuming part of the context structure associated to the GKB algorithm.
[in,out] | ctx | pointer to the context structure to clean |
void cs_saddle_solver_context_gkb_create | ( | cs_saddle_solver_t * | solver | ) |
Create and initialize the context structure for the GKB algorithm.
[in,out] | solver | pointer to a saddle-point solver structure |
void cs_saddle_solver_context_gkb_free | ( | cs_saddle_solver_context_gkb_t ** | p_ctx | ) |
Free the context structure associated to the GKB algorithm.
[in,out] | p_ctx | double pointer to the context structure to free |
void cs_saddle_solver_context_notay_create | ( | cs_saddle_solver_t * | solver | ) |
Create and initialize the context structure for the algorithm relying on the Notay's algebraic transformation.
[in,out] | solver | pointer to a saddle-point solver structure |
void cs_saddle_solver_context_simple_clean | ( | cs_saddle_solver_context_simple_t * | ctx | ) |
Free main memory consuming part of the context structure associated to a SIMPLE algorithm.
[in,out] | ctx | pointer to the context structure to clean |
void cs_saddle_solver_context_simple_create | ( | cs_lnum_t | b22_max_size, |
cs_saddle_solver_t * | solver | ||
) |
Create and initialize the context structure for an algorithm related to the SIMPLE algorithm.
[in] | b22_max_size | max. size for the second part of unknows |
[in,out] | solver | pointer to a saddle-point solver structure |
Create and initialize the context structure for an algorithm related to the SIMPLE algorithm.
[in] | b22_max_size | max. size for the second part of unknows |
[in,out] | solver | pointer to a saddle-point solver structure |
void cs_saddle_solver_context_simple_free | ( | cs_saddle_solver_context_simple_t ** | p_ctx | ) |
Free the context structure associated to a Uzawa-CG algorithm.
[in,out] | p_ctx | double pointer to the context structure to free |
Free the context structure associated to a Uzawa-CG algorithm.
[in,out] | p_ctx | double pointer to the context structure to free |
void cs_saddle_solver_context_uzawa_cg_clean | ( | cs_saddle_solver_context_uzawa_cg_t * | ctx | ) |
Free main memory consuming part of the context structure associated to a Uzawa-CG algorithm.
[in,out] | ctx | pointer to the context structure to clean |
void cs_saddle_solver_context_uzawa_cg_create | ( | cs_lnum_t | b22_max_size, |
cs_saddle_solver_t * | solver | ||
) |
Create and initialize the context structure for an algorithm related to the Uzawa-CG algorithm.
[in] | b22_max_size | max. size for the second part of unknows |
[in,out] | solver | pointer to a saddle-point solver structure |
void cs_saddle_solver_context_uzawa_cg_free | ( | cs_saddle_solver_context_uzawa_cg_t ** | p_ctx | ) |
Free the context structure associated to a Uzawa-CG algorithm.
[in,out] | p_ctx | double pointer to the context structure to free |
void cs_saddle_solver_finalize | ( | void | ) |
Free all remaining structures related to saddle-point solvers.
void cs_saddle_solver_free | ( | cs_saddle_solver_t ** | p_solver | ) |
Free a cs_saddle_solver_t structure.
[in,out] | p_solver | double pointer to the structure to free |
void cs_saddle_solver_gcr | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | x1, | ||
cs_real_t * | x2 | ||
) |
Apply the GCR algorithm to a saddle point problem (the system is stored in a hybrid way). The stride is equal to 1 for the matrix in the (1,1) block (db_size[3] = 1).
This algorithm is taken from 2010 Notay's paper: "An aggregation-based algebraic multigrid method" ETNA (vol. 37)
[in,out] | solver | pointer to a saddle_point solver structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
int cs_saddle_solver_get_n_systems | ( | void | ) |
Retrieve the number of saddle-point systems which have been added.
void cs_saddle_solver_gkb_inhouse | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | x1, | ||
cs_real_t * | x2 | ||
) |
Apply the GKB algorithm to a saddle point problem (the system is stored in a hybrid way). The stride is equal to 1 for the matrix in the (1,1) block (db_size[3] = 1).
[in,out] | solver | pointer to a saddle_point solver structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
void cs_saddle_solver_log_monitoring | ( | void | ) |
Log the monitoring performance of all saddle-point systems.
Log the monitoring performance of all saddle-point systems.
cs_real_t * cs_saddle_solver_m11_inv_lumped | ( | cs_saddle_solver_t * | solver, |
const cs_matrix_t * | m11, | ||
const cs_range_set_t * | b11_rset, | ||
cs_sles_t * | xtra_sles, | ||
int * | n_iter | ||
) |
Retrieve the lumped matrix 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 scatter view.
[in] | solver | solver for saddle-point problems |
[in] | m11 | matrix related to the (1,1) block |
[in] | b11_rset | range set structure for the (1,1) block |
[in,out] | xtra_sles | pointer to an extra SLES structure |
[out] | n_iter | number of iterations for this operation |
void cs_saddle_solver_m12_multiply_scalar | ( | cs_lnum_t | n2_elts, |
const cs_real_t * | x2, | ||
const cs_adjacency_t * | m21_adj, | ||
const cs_real_t * | m21_val, | ||
cs_real_t * | m12x2 | ||
) |
Compute the resulting vector of the operation m12*x2 The stride is equal to 1 for the operator m21 (unassembled) x2 corresponds to a "scatter" view.
This is an update operation. Be careful that the resulting array has to be initialized.
[in] | n2_elts | number of elements for x2 (not DoFs) |
[in] | x2 | array for the second set |
[in] | m21_adj | adjacency related to the M21 operator |
[in] | m21_val | array associated to the M21 operator (unassembled) |
[in,out] | m12x2 | resulting array (have to be allocated) |
void cs_saddle_solver_m12_multiply_vector | ( | cs_lnum_t | n2_dofs, |
const cs_real_t * | x2, | ||
const cs_adjacency_t * | m21_adj, | ||
const cs_real_t * | m21_val, | ||
cs_real_t * | m12x2 | ||
) |
Compute the resulting vector of the operation m12*x2 The stride is equal to 3 for the operator m21 (unassembled) x2 corresponds to a "scatter" view.
This is an update operation. Be careful that the resulting array has to be initialized.
[in] | n2_dofs | number of DoFs for x2 |
[in] | x2 | array for the second set |
[in] | m21_adj | adjacency related to the M21 operator |
[in] | m21_val | array associated to the M21 operator (unassembled) |
[in,out] | m12x2 | resulting array (have to be allocated) |
void cs_saddle_solver_m21_multiply_scalar | ( | cs_lnum_t | n2_dofs, |
const cs_real_t * | x1, | ||
const cs_adjacency_t * | m21_adj, | ||
const cs_real_t * | m21_val, | ||
cs_real_t * | m21x1 | ||
) |
Compute the resulting vector of the operation m21*x1 The stride is equal to 1 for the operator m21 operator x1 corresponds to a "scatter" view.
[in] | n2_dofs | number of (scatter) DoFs for (2,2)-block |
[in] | x1 | array for the first part |
[in] | m21_adj | adjacency related to the M21 operator |
[in] | m21_val | values associated to the M21 operator (unassembled) |
[in,out] | m21x1 | resulting vector (have to be allocated) |
void cs_saddle_solver_m21_multiply_vector | ( | cs_lnum_t | n2_dofs, |
const cs_real_t * | x1, | ||
const cs_adjacency_t * | m21_adj, | ||
const cs_real_t * | m21_val, | ||
cs_real_t * | m21x1 | ||
) |
Compute the resulting vector of the operation m21*x1 The stride is equal to 3 for the operator m21 operator x1 corresponds to a "scatter" view.
[in] | n2_dofs | number of (scatter) DoFs for (2,2)-block |
[in] | x1 | array for the first part |
[in] | m21_adj | adjacency related to the M21 operator |
[in] | m21_val | values associated to the M21 operator (unassembled) |
[in,out] | m21x1 | resulting vector (have to be allocated) |
void cs_saddle_solver_minres | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | x1, | ||
cs_real_t * | x2 | ||
) |
Apply the MINRES algorithm to a saddle point problem (the system is stored in a hybrid way). The stride is equal to 1 for the matrix (db_size[3] = 1) and the vector.
[in,out] | solver | pointer to a saddle_point solver structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
void cs_saddle_solver_notay | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | x1, | ||
cs_real_t * | x2 | ||
) |
Apply the Notay's transformation algorithm to solve a saddle point problem (the system is stored in a monolithic way).
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
void cs_saddle_solver_simple | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | x1, | ||
cs_real_t * | x2 | ||
) |
Apply the SIMPLE-like algorithm to solve a saddle point problem.
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
void cs_saddle_solver_sles_full_system | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | x1, | ||
cs_real_t * | x2 | ||
) |
Apply an (external) solver to solve a saddle point problem (the system is stored in a monolithic way).
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |
void cs_saddle_solver_update_monitoring | ( | cs_saddle_solver_t * | solver, |
unsigned | n_iter | ||
) |
Update the current monitoring state with n_iter.
[in,out] | solver | pointer to a saddle solver structure |
[in] | n_iter | number of iterations needed for a new call |
void cs_saddle_solver_uzawa_cg | ( | cs_saddle_solver_t * | solver, |
cs_real_t * | x1, | ||
cs_real_t * | x2 | ||
) |
Apply the Uzawa-CG algorithm to solve a saddle point problem (the system is stored in a hybrid way). This algorithm is based on Koko's paper "Uzawa conjugate gradient method for the Stokes problem: Matlab implementation with P1-iso-P2/ P1 finite element".
[in,out] | solver | pointer to a cs_saddle_solver_t structure |
[in,out] | x1 | array for the first part |
[in,out] | x2 | array for the second part |