8.3
general documentation
cs_saddle_solver.h File Reference
#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"
+ Include dependency graph for cs_saddle_solver.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_tcs_saddle_solver_by_id (int id)
 Get a pointer to saddle-point solver from its id. More...
 
cs_saddle_solver_tcs_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_tcs_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 Documentation

◆ cs_saddle_solver_matvec_t

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.

Parameters
[in]n2_dofsnumber of (scatter) DoFs for the (2,2)-block
[in]vecarray of values
[in]mat_adjadjacency related to the matrix operator
[in]mat_valvalues associated to the matrix operator
[in,out]matvecresulting vector (have to be allocated)

Function Documentation

◆ cs_saddle_solver_add()

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.

Parameters
[in]n1_eltsnumber of elements associated to the (1,1)-block
[in]n1_dofs_by_eltnumber of DoFs by elements in the (1,1)-block
[in]n2_eltsnumber of elements associated to the (2,2)-block
[in]n2_dofs_by_eltnumber of DoFs by elements in the (2,2)-block
[in]saddlepset of parameters for the saddle-point solver
[in]shpointer to a system helper structure
[in]main_slespointer to the main SLES structure related to this saddle-point problem
Returns
a pointer to the new allocated structure

◆ cs_saddle_solver_alu_incr()

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.

Parameters
[in,out]solverpointer to a saddle_point solver structure
[in,out]x1array for the first part
[in,out]x2array for the second part

◆ cs_saddle_solver_by_id()

cs_saddle_solver_t * cs_saddle_solver_by_id ( int  id)

Get a pointer to saddle-point solver from its id.

Parameters
[in]idid of the saddle-point system
Returns
a pointer to a saddle-point solver structure

Get a pointer to saddle-point solver from its id.

Parameters
[in]idid of the saddle-point system
Returns
a pointer to a saddle-point solver structure

◆ cs_saddle_solver_clean()

void cs_saddle_solver_clean ( cs_saddle_solver_t solver)

Free/reset only a part of a cs_saddle_solver_t structure.

Parameters
[in,out]solverdouble pointer to the structure to free

◆ cs_saddle_solver_context_alu_clean()

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.

Parameters
[in,out]ctxpointer to the context structure to clean

◆ cs_saddle_solver_context_alu_create()

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.

Parameters
[in,out]solverpointer to a saddle-point solver structure

◆ cs_saddle_solver_context_alu_free()

void cs_saddle_solver_context_alu_free ( cs_saddle_solver_context_alu_t **  p_ctx)

Free the context structure associated to an ALU algorithm.

Parameters
[in,out]p_ctxdouble pointer to the context structure to free

◆ cs_saddle_solver_context_block_pcd_clean()

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.

Parameters
[in,out]ctxpointer 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.

Parameters
[in,out]ctxpointer to the context structure to free

◆ cs_saddle_solver_context_block_pcd_create()

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.

Parameters
[in]b22_max_sizemax. size for the second part of unknows
[in,out]solverpointer to a saddle-point solver structure

◆ cs_saddle_solver_context_block_pcd_free()

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.

Parameters
[in,out]p_ctxdouble pointer to the context structure to free

◆ cs_saddle_solver_context_gkb_clean()

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.

Parameters
[in,out]ctxpointer to the context structure to clean

◆ cs_saddle_solver_context_gkb_create()

void cs_saddle_solver_context_gkb_create ( cs_saddle_solver_t solver)

Create and initialize the context structure for the GKB algorithm.

Parameters
[in,out]solverpointer to a saddle-point solver structure

◆ cs_saddle_solver_context_gkb_free()

void cs_saddle_solver_context_gkb_free ( cs_saddle_solver_context_gkb_t **  p_ctx)

Free the context structure associated to the GKB algorithm.

Parameters
[in,out]p_ctxdouble pointer to the context structure to free

◆ cs_saddle_solver_context_notay_create()

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.

Parameters
[in,out]solverpointer to a saddle-point solver structure

◆ cs_saddle_solver_context_simple_clean()

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.

Parameters
[in,out]ctxpointer to the context structure to clean

◆ cs_saddle_solver_context_simple_create()

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.

Parameters
[in]b22_max_sizemax. size for the second part of unknows
[in,out]solverpointer to a saddle-point solver structure

Create and initialize the context structure for an algorithm related to the SIMPLE algorithm.

Parameters
[in]b22_max_sizemax. size for the second part of unknows
[in,out]solverpointer to a saddle-point solver structure

◆ cs_saddle_solver_context_simple_free()

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.

Parameters
[in,out]p_ctxdouble pointer to the context structure to free

Free the context structure associated to a Uzawa-CG algorithm.

Parameters
[in,out]p_ctxdouble pointer to the context structure to free

◆ cs_saddle_solver_context_uzawa_cg_clean()

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.

Parameters
[in,out]ctxpointer to the context structure to clean

◆ cs_saddle_solver_context_uzawa_cg_create()

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.

Parameters
[in]b22_max_sizemax. size for the second part of unknows
[in,out]solverpointer to a saddle-point solver structure

◆ cs_saddle_solver_context_uzawa_cg_free()

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.

Parameters
[in,out]p_ctxdouble pointer to the context structure to free

◆ cs_saddle_solver_finalize()

void cs_saddle_solver_finalize ( void  )

Free all remaining structures related to saddle-point solvers.

◆ cs_saddle_solver_free()

void cs_saddle_solver_free ( cs_saddle_solver_t **  p_solver)

Free a cs_saddle_solver_t structure.

Parameters
[in,out]p_solverdouble pointer to the structure to free

◆ cs_saddle_solver_gcr()

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)

Parameters
[in,out]solverpointer to a saddle_point solver structure
[in,out]x1array for the first part
[in,out]x2array for the second part

◆ cs_saddle_solver_get_n_systems()

int cs_saddle_solver_get_n_systems ( void  )

Retrieve the number of saddle-point systems which have been added.

Returns
the current number of systems

◆ cs_saddle_solver_gkb_inhouse()

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

Parameters
[in,out]solverpointer to a saddle_point solver structure
[in,out]x1array for the first part
[in,out]x2array for the second part

◆ cs_saddle_solver_log_monitoring()

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_saddle_solver_m11_inv_lumped()

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.

Parameters
[in]solversolver for saddle-point problems
[in]m11matrix related to the (1,1) block
[in]b11_rsetrange set structure for the (1,1) block
[in,out]xtra_slespointer to an extra SLES structure
[out]n_iternumber of iterations for this operation
Returns
a pointer to the computed array (scatter view)

◆ cs_saddle_solver_m12_multiply_scalar()

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.

Parameters
[in]n2_eltsnumber of elements for x2 (not DoFs)
[in]x2array for the second set
[in]m21_adjadjacency related to the M21 operator
[in]m21_valarray associated to the M21 operator (unassembled)
[in,out]m12x2resulting array (have to be allocated)

◆ cs_saddle_solver_m12_multiply_vector()

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.

Parameters
[in]n2_dofsnumber of DoFs for x2
[in]x2array for the second set
[in]m21_adjadjacency related to the M21 operator
[in]m21_valarray associated to the M21 operator (unassembled)
[in,out]m12x2resulting array (have to be allocated)

◆ cs_saddle_solver_m21_multiply_scalar()

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.

Parameters
[in]n2_dofsnumber of (scatter) DoFs for (2,2)-block
[in]x1array for the first part
[in]m21_adjadjacency related to the M21 operator
[in]m21_valvalues associated to the M21 operator (unassembled)
[in,out]m21x1resulting vector (have to be allocated)

◆ cs_saddle_solver_m21_multiply_vector()

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.

Parameters
[in]n2_dofsnumber of (scatter) DoFs for (2,2)-block
[in]x1array for the first part
[in]m21_adjadjacency related to the M21 operator
[in]m21_valvalues associated to the M21 operator (unassembled)
[in,out]m21x1resulting vector (have to be allocated)

◆ cs_saddle_solver_minres()

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.

Parameters
[in,out]solverpointer to a saddle_point solver structure
[in,out]x1array for the first part
[in,out]x2array for the second part

◆ cs_saddle_solver_notay()

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

Parameters
[in,out]solverpointer to a cs_saddle_solver_t structure
[in,out]x1array for the first part
[in,out]x2array for the second part

◆ cs_saddle_solver_simple()

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.

Parameters
[in,out]solverpointer to a cs_saddle_solver_t structure
[in,out]x1array for the first part
[in,out]x2array for the second part

◆ cs_saddle_solver_sles_full_system()

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

Parameters
[in,out]solverpointer to a cs_saddle_solver_t structure
[in,out]x1array for the first part
[in,out]x2array for the second part

◆ cs_saddle_solver_update_monitoring()

void cs_saddle_solver_update_monitoring ( cs_saddle_solver_t solver,
unsigned  n_iter 
)

Update the current monitoring state with n_iter.

Parameters
[in,out]solverpointer to a saddle solver structure
[in]n_iternumber of iterations needed for a new call

◆ cs_saddle_solver_uzawa_cg()

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".

Parameters
[in,out]solverpointer to a cs_saddle_solver_t structure
[in,out]x1array for the first part
[in,out]x2array for the second part