8.3
general documentation
cs_param_saddle.h File Reference

Handle the settings of saddle-point systems. These systems arise from the monolithic coupling of the Navier-Stokes equations or in mixed formulation of scalar-valued equations. More...

#include "cs_param_sles.h"
+ Include dependency graph for cs_param_saddle.h:

Go to the source code of this file.

Data Structures

struct  cs_param_saddle_t
 Structure storing all metadata related to the resolution of a saddle-point linear system. A saddle-point system is depicted as. More...
 
struct  cs_param_saddle_context_alu_t
 
struct  cs_param_saddle_context_block_krylov_t
 
struct  cs_param_saddle_context_gkb_t
 
struct  cs_param_saddle_context_notay_t
 
struct  cs_param_saddle_context_uzacg_t
 
struct  cs_param_saddle_context_simple_t
 

Enumerations

enum  cs_param_saddle_precond_t {
  CS_PARAM_SADDLE_PRECOND_NONE , CS_PARAM_SADDLE_PRECOND_DIAG , CS_PARAM_SADDLE_PRECOND_LOWER , CS_PARAM_SADDLE_PRECOND_SGS ,
  CS_PARAM_SADDLE_PRECOND_UPPER , CS_PARAM_SADDLE_PRECOND_UZAWA , CS_PARAM_SADDLE_N_PRECOND
}
 Type of preconditioner used to solve a saddle-point system. Up to now, this happens only in two cases: (1) with CDO cell-based schemes and (2) with the Stokes or Navier-Stokes equations with CDO face-based schemes and a monolithic approach (fully coupled) velocity-pressure coupling. More...
 
enum  cs_param_saddle_solver_t {
  CS_PARAM_SADDLE_SOLVER_NONE , CS_PARAM_SADDLE_SOLVER_ALU , CS_PARAM_SADDLE_SOLVER_FGMRES , CS_PARAM_SADDLE_SOLVER_GCR ,
  CS_PARAM_SADDLE_SOLVER_GKB , CS_PARAM_SADDLE_SOLVER_MINRES , CS_PARAM_SADDLE_SOLVER_MUMPS , CS_PARAM_SADDLE_SOLVER_NOTAY_TRANSFORM ,
  CS_PARAM_SADDLE_SOLVER_SIMPLE , CS_PARAM_SADDLE_SOLVER_UZAWA_CG , CS_PARAM_SADDLE_N_SOLVERS
}
 Type of solver used to solve a saddle-point system. Up to now, this happens only with CDO cell-based schemes or when solving the fully coupled Navier-Stokes system (monolithic) with CDO face-based schemes. A saddle-point system is an indefinite system. More...
 
enum  cs_param_saddle_schur_approx_t {
  CS_PARAM_SADDLE_SCHUR_NONE , CS_PARAM_SADDLE_SCHUR_DIAG_INVERSE , CS_PARAM_SADDLE_SCHUR_IDENTITY , CS_PARAM_SADDLE_SCHUR_LUMPED_INVERSE ,
  CS_PARAM_SADDLE_SCHUR_MASS_SCALED , CS_PARAM_SADDLE_SCHUR_MASS_SCALED_DIAG_INVERSE , CS_PARAM_SADDLE_SCHUR_MASS_SCALED_LUMPED_INVERSE , CS_PARAM_SADDLE_N_SCHUR_APPROX
}
 Strategy to build the Schur complement approximation. This appears in block preconditioning or Uzawa algorithms when a fully coupled (also called monolithic) approach) is used. More...
 

Functions

void cs_param_saddle_set_restart_range (cs_param_saddle_t *saddlep, int restart_range)
 Set the number of iterations to store before starting a Krylov solver. More...
 
void cs_param_saddle_set_notay_scaling (cs_param_saddle_t *saddlep, double scaling_coef)
 Set the scaling coefficient used in the Notay's transformation devised in "Algebraic multigrid for Stokes equations" SIAM J. Sci. Comput. Vol. 39 (5), 2017 In this article, this scaling is denoted by alpha. More...
 
void cs_param_saddle_set_augmentation_coef (cs_param_saddle_t *saddlep, double coef)
 Set the scaling in front of the augmentation term when an ALU, a GKB or a block Krylov algorithm is considered. More...
 
double cs_param_saddle_get_augmentation_coef (const cs_param_saddle_t *saddlep)
 Get the scaling coefficient in front of the augmentation term when an ALU, GKB or block Krylov algorithm is considered. More...
 
const char * cs_param_saddle_get_type_name (cs_param_saddle_solver_t type)
 Retrieve the name of the type of saddle-point solver. More...
 
cs_param_saddle_tcs_param_saddle_create (void)
 Create a cs_param_saddle_t structure No solver is set by default. More...
 
void cs_param_saddle_free (cs_param_saddle_t **p_saddlep)
 Free the structure storing the parameter settings for a saddle-point system. More...
 
const char * cs_param_saddle_get_name (const cs_param_saddle_t *saddlep)
 Retrieve the name of the saddle-point solver. More...
 
void cs_param_saddle_set_name (const char *name, cs_param_saddle_t *saddlep)
 Set the name of the saddle-point system. More...
 
void cs_param_saddle_set_block11_sles_param (cs_param_saddle_t *saddlep, const cs_param_sles_t *block11_slesp)
 Assign the cs_param_sles_t structure (shared) related to the (1,1)-block to the structure managing the resolution of the saddle-point problems. More...
 
int cs_param_saddle_set_precond (const char *keyval, cs_param_saddle_t *saddlep)
 Set the type of preconditioning to apply for this saddle-point system. More...
 
int cs_param_saddle_set_schur_approx (const char *keyval, cs_param_saddle_t *saddlep)
 Set the type of Schur approximation to apply to this saddle-point system. More...
 
int cs_param_saddle_set_solver_class (const char *keyval, cs_param_saddle_t *saddlep)
 Set the class of solver to apply for this saddle-point system. More...
 
int cs_param_saddle_set_solver (const char *keyval, cs_param_saddle_t *saddlep)
 Set the type of solver to apply for this saddle-point system. More...
 
void cs_param_saddle_try_init_schur_sles_param (cs_param_saddle_t *saddlep)
 Initialize a cs_param_sles_t structure for the Schur approximation nested inside a cs_param_saddle_t structure. By default, this member is not allocated. Do nothing if the related structure is already allocated. More...
 
cs_param_sles_tcs_param_saddle_get_schur_sles_param (const cs_param_saddle_t *saddlep)
 Get the pointer to the set of parameters to handle a SLES. This SLES is associated to the approximation of the Schur complement. This is only useful for solving a saddle-point problem relying on an elaborated approximation of the Schur complement. More...
 
cs_param_sles_tcs_param_saddle_get_xtra_sles_param (const cs_param_saddle_t *saddlep)
 Get the pointer to the set of parameters to handle a SLES. This SLES is associated to an extra-operation specific to a saddle-point solver It returns a non-null pointer only for some sadlle-point solver relying on a more elaborated Schur complement approximation. More...
 
cs_param_sles_tcs_param_saddle_get_init_sles_param (const cs_param_saddle_t *saddlep)
 Get the pointer to the set of parameters to handle a SLES. This SLES is associated to the initial saddle-point problem. It returns a non null pointer only for some sadlle-point solver. More...
 
void cs_param_saddle_copy (const cs_param_saddle_t *ref, cs_param_saddle_t *dest)
 Copy a cs_param_saddle_t structure from ref to dest. More...
 
void cs_param_saddle_log (const cs_param_saddle_t *saddlep)
 Log the setup information for the given cs_param_saddle_t structure. More...
 

Detailed Description

Handle the settings of saddle-point systems. These systems arise from the monolithic coupling of the Navier-Stokes equations or in mixed formulation of scalar-valued equations.

Enumeration Type Documentation

◆ cs_param_saddle_precond_t

Type of preconditioner used to solve a saddle-point system. Up to now, this happens only in two cases: (1) with CDO cell-based schemes and (2) with the Stokes or Navier-Stokes equations with CDO face-based schemes and a monolithic approach (fully coupled) velocity-pressure coupling.

Enumerator
CS_PARAM_SADDLE_PRECOND_NONE 

No preconditioner to apply.

CS_PARAM_SADDLE_PRECOND_DIAG 

A block-diagonal preconditioner is used. The (1,1)-block is an approximation of the (1,1)-block in the saddle-point system based on a cheap resolution. The parameter settings for this resolution relies on the structure cs_param_sles_t. The (2,2)-block is approximated by the Schur approximation given (by default the identity).

CS_PARAM_SADDLE_PRECOND_LOWER 

A 2x2 block matrix is used as preconditioner with the (1,2)-block fills with zero. The (1,1)-block is an approximation of the (1,1)-block in the saddle-point system based on a cheap resolution. The parameter settings for this resolution relies on the structure cs_param_sles_t. The (2,2)-block is approximated by the Schur approximation given (by default the identity).

CS_PARAM_SADDLE_PRECOND_SGS 

A symmetric Gauss-Seidel 2x2 block preconditioner is used. The (1,1)-block is an approximation of the (1,1)-block in the saddle-point system based on a cheap resolution. The parameter settings for this resolution relies on the structure cs_param_sles_t The (2,2)-block is approximated by the Schur approximation given (by default the identity).

CS_PARAM_SADDLE_PRECOND_UPPER 

A 2x2 block matrix is used as preconditioner with the (2,1)-block fills with zero. The (1,1)-block is an approximation of the (1,1)-block in the saddle-point system based on a cheap resolution. The parameter settings for this resolution relies on the structure cs_param_sles_t The (2,2)-block is approximated by the Schur approximation given (by default the identity).

CS_PARAM_SADDLE_PRECOND_UZAWA 

An Uzawa-like 2x2 block preconditioner. One needs a Schur approximation.

CS_PARAM_SADDLE_N_PRECOND 

◆ cs_param_saddle_schur_approx_t

Strategy to build the Schur complement approximation. This appears in block preconditioning or Uzawa algorithms when a fully coupled (also called monolithic) approach) is used.

\[ \begin{bmatrix} A&B^t\\ B&O  \end{bmatrix}\]

The exact Schur complement is then

\[ S = -B \cdot A^{-1} \cdot B \]

Enumerator
CS_PARAM_SADDLE_SCHUR_NONE 

There is no Schur complement approximation.

CS_PARAM_SADDLE_SCHUR_DIAG_INVERSE 

The Schur complement approximation is defined as

\[ S \approx -B \cdot diag(A)^{-1} \cdot B^t \]

CS_PARAM_SADDLE_SCHUR_IDENTITY 

The Schur complement approximation is simply the identity matrix

CS_PARAM_SADDLE_SCHUR_LUMPED_INVERSE 

The Schur complement approximation is defined as

\[ B \cdot lumped(A^{-1}) \cdot B^t \]

where $x=lumped(A^{-1})$ results from $A.x = \bf{1}$ ( $\bf{1}$ is the array fills with 1 in each entry)

CS_PARAM_SADDLE_SCHUR_MASS_SCALED 

The Schur complement approximation is simply a scaled diagonal mass matrix related to the (2,2)-block

CS_PARAM_SADDLE_SCHUR_MASS_SCALED_DIAG_INVERSE 

The Schur complement approximation is defined as

\[ S \approx \alpha M_{22} + \frac{1}{dt} B.diag(A)^{-1}.B^t \]

where $M_{22}$ is the mass matrix related to the (2,2)-block

CS_PARAM_SADDLE_SCHUR_MASS_SCALED_LUMPED_INVERSE 

The Schur complement approximation is defined as

\[ S \approx \alpha \cdot M_{22} + \frac{1}{dt} B\cdot lumped(A^{-1})
 \cdot B^t \]

where $ M_{22} $ is the mass matrix related to the (2,2) block and where $x=lumped(A^{-1})$ results from $A.x = \bf{1}$ ( $\bf{1}$ is the array fills with 1 in each entry)

CS_PARAM_SADDLE_N_SCHUR_APPROX 

◆ cs_param_saddle_solver_t

Type of solver used to solve a saddle-point system. Up to now, this happens only with CDO cell-based schemes or when solving the fully coupled Navier-Stokes system (monolithic) with CDO face-based schemes. A saddle-point system is an indefinite system.

Enumerator
CS_PARAM_SADDLE_SOLVER_NONE 

No solver defined. No saddle-point to solve.

CS_PARAM_SADDLE_SOLVER_ALU 

In-house solver based on the Uzawa algorithm with an Augmented Lagrangian acceleration. In this case, one has to specify the scaling coefficient for the augmented system (in the current case, a grad-div operator)

CS_PARAM_SADDLE_SOLVER_FGMRES 

Flexible variant of the GMRES Krylov iterative solver. This solver can be applied to a general indefinite systems. Up to now, this solver can be called only through the PETSc library.

CS_PARAM_SADDLE_SOLVER_GCR 

GCR = Genreralized Conjugate Residual. This iterative solver can be applied to a general indefinite systems. It relies on a specific storage of the saddle-point system: (1,1)-block is assembled and the (2,1)-block is unassembled. The (1,2) is not stored since this is the transposition of the (2,1)-block. This solver is an in-house version of the Notay's variant for the GCR. One can specify the solver for the (1,1)-block preconditioner and also the type of Schur approximation for the (2,2)-block.

CS_PARAM_SADDLE_SOLVER_GKB 

GKB = Golub-Kahan bi-diagonalization. This solver can be applied to a symmetrice saddle-point system (or very slightly unsymmetric system). The stopping crietrion is different from the other solvers since it relies on an approximation of the residual in the energy norm (evaluation on-the-fly). This solver is available throught PETSc or directly with code_saturne. Moreover, one can add an augmentation term as in the ALU algorithm.

CS_PARAM_SADDLE_SOLVER_MINRES 

MinRes = Minimal Residual algorithm. This iterative solver can be applied to an indefinite symmetric systems. This is like a CG solver for symmetric saddle-point systems. This solver relies on a specific storage of the saddle-point system: (1,1)-block is assembled and the (2,1)-block is unassembled. The (1,2) is not stored since this the transposition of the (2,1)-block. One can specify the solver for the (1,1)-block preconditioner and also the type of Schur approximation for the (2,2)-block.

CS_PARAM_SADDLE_SOLVER_MUMPS 

Sparse direct solver based on the external library called MUMPS. This iterative solver can be applied to a general indefinite systems. Be careful that according to the boundary conditions. The pressure can be defined up to a constant. In this situation, MUMPS can breakdown since there is a non empty kernel.

CS_PARAM_SADDLE_SOLVER_NOTAY_TRANSFORM 

Transform the saddle-point problem into an equivalent system without a zero block for the (2,2) block. This transformation allows one to consider standard preconditionner and iterative Krylov solver. This is an experimental feature (only tested in sequential run)

CS_PARAM_SADDLE_SOLVER_SIMPLE 
CS_PARAM_SADDLE_SOLVER_UZAWA_CG 

Resolution using an uzawa algorithm optimized using a conjugate gradient reformulation. Two systems are solved at each iteration (one related to the velocity block, one related to the Schur complement approximation - size of the pressure space). The first one is specified using the cs_sles_param_t structure related to the (1,1)-block and the second system is specified using the cs_sles_param_t structure related to the Schur approximation.

CS_PARAM_SADDLE_N_SOLVERS 

Function Documentation

◆ cs_param_saddle_copy()

void cs_param_saddle_copy ( const cs_param_saddle_t ref,
cs_param_saddle_t dest 
)

Copy a cs_param_saddle_t structure from ref to dest.

Parameters
[in]refreference structure to be copied
[in,out]destdestination structure

◆ cs_param_saddle_create()

cs_param_saddle_t * cs_param_saddle_create ( void  )

Create a cs_param_saddle_t structure No solver is set by default.

Returns
a pointer to the new cs_param_saddle_t structure

◆ cs_param_saddle_free()

void cs_param_saddle_free ( cs_param_saddle_t **  p_saddlep)

Free the structure storing the parameter settings for a saddle-point system.

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

◆ cs_param_saddle_get_augmentation_coef()

double cs_param_saddle_get_augmentation_coef ( const cs_param_saddle_t saddlep)

Get the scaling coefficient in front of the augmentation term when an ALU, GKB or block Krylov algorithm is considered.

Parameters
[in]saddlepset of parameters for solving a saddle-point
Returns
0 if not relevant or the value of the augmentation coefficient

Get the scaling coefficient in front of the augmentation term when an ALU, GKB or block Krylov algorithm is considered.

Parameters
[in]saddlepset of parameters for solving a saddle-point
Returns
0 if not relevant or the value of the augmentation coefficient

◆ cs_param_saddle_get_init_sles_param()

cs_param_sles_t * cs_param_saddle_get_init_sles_param ( const cs_param_saddle_t saddlep)

Get the pointer to the set of parameters to handle a SLES. This SLES is associated to the initial saddle-point problem. It returns a non null pointer only for some sadlle-point solver.

Parameters
[in]saddleppointer to a cs_param_saddle_t structure
Returns
a pointer to a cs_param_sles_t structure

◆ cs_param_saddle_get_name()

const char * cs_param_saddle_get_name ( const cs_param_saddle_t saddlep)

Retrieve the name of the saddle-point solver.

Parameters
[in]saddleppointer to a set of saddle-point parameters
Returns
a string

◆ cs_param_saddle_get_schur_sles_param()

cs_param_sles_t * cs_param_saddle_get_schur_sles_param ( const cs_param_saddle_t saddlep)

Get the pointer to the set of parameters to handle a SLES. This SLES is associated to the approximation of the Schur complement. This is only useful for solving a saddle-point problem relying on an elaborated approximation of the Schur complement.

Parameters
[in]saddleppointer to a cs_param_saddle_t structure
Returns
a pointer to a cs_param_sles_t structure

◆ cs_param_saddle_get_type_name()

const char * cs_param_saddle_get_type_name ( cs_param_saddle_solver_t  type)

Retrieve the name of the type of saddle-point solver.

Parameters
[in]typetype of saddle-point solver
Returns
a string

◆ cs_param_saddle_get_xtra_sles_param()

cs_param_sles_t * cs_param_saddle_get_xtra_sles_param ( const cs_param_saddle_t saddlep)

Get the pointer to the set of parameters to handle a SLES. This SLES is associated to an extra-operation specific to a saddle-point solver It returns a non-null pointer only for some sadlle-point solver relying on a more elaborated Schur complement approximation.

Parameters
[in]saddleppointer to a cs_param_saddle_t structure
Returns
a pointer to a cs_param_sles_t structure

◆ cs_param_saddle_log()

void cs_param_saddle_log ( const cs_param_saddle_t saddlep)

Log the setup information for the given cs_param_saddle_t structure.

Parameters
[in]saddleppointer to the structure

◆ cs_param_saddle_set_augmentation_coef()

void cs_param_saddle_set_augmentation_coef ( cs_param_saddle_t saddlep,
double  coef 
)

Set the scaling in front of the augmentation term when an ALU, a GKB or a block Krylov algorithm is considered.

Parameters
[in,out]saddlepset of parameters for solving a saddle-point
[in]coefvalue of the scaling coefficient

Set the scaling in front of the augmentation term when an ALU, a GKB or a block Krylov algorithm is considered.

Parameters
[in,out]saddlepset of parameters for solving a saddle-point
[in]coefvalue of the scaling coefficient

◆ cs_param_saddle_set_block11_sles_param()

void cs_param_saddle_set_block11_sles_param ( cs_param_saddle_t saddlep,
const cs_param_sles_t block11_slesp 
)

Assign the cs_param_sles_t structure (shared) related to the (1,1)-block to the structure managing the resolution of the saddle-point problems.

Parameters
[in,out]saddleppointer to the structure to update
[in]block11_slespset of parameters for the (1,1) block

◆ cs_param_saddle_set_name()

void cs_param_saddle_set_name ( const char *  name,
cs_param_saddle_t saddlep 
)

Set the name of the saddle-point system.

Parameters
[in]namename associated to this saddle-point system
[in,out]saddleppointer to the structure to update

◆ cs_param_saddle_set_notay_scaling()

void cs_param_saddle_set_notay_scaling ( cs_param_saddle_t saddlep,
double  scaling_coef 
)

Set the scaling coefficient used in the Notay's transformation devised in "Algebraic multigrid for Stokes equations" SIAM J. Sci. Comput. Vol. 39 (5), 2017 In this article, this scaling is denoted by alpha.

Parameters
[in,out]saddlepset of parameters for solving a saddle-point
[in]scaling_coefvalue of the scaling coefficient

◆ cs_param_saddle_set_precond()

int cs_param_saddle_set_precond ( const char *  keyval,
cs_param_saddle_t saddlep 
)

Set the type of preconditioning to apply for this saddle-point system.

Parameters
[in]keyvalvalue of the key for the preconditioner
[in,out]saddleppointer to the structure to update
Returns
0 if no error detected, > 0 otherwise

◆ cs_param_saddle_set_restart_range()

void cs_param_saddle_set_restart_range ( cs_param_saddle_t saddlep,
int  restart_range 
)

Set the number of iterations to store before starting a Krylov solver.

Parameters
[in,out]saddlepset of parameters for solving a saddle-point
[in]restart_rangenumber of directions

◆ cs_param_saddle_set_schur_approx()

int cs_param_saddle_set_schur_approx ( const char *  keyval,
cs_param_saddle_t saddlep 
)

Set the type of Schur approximation to apply to this saddle-point system.

Parameters
[in]keyvalvalue of the key for the schur approx.
[in,out]saddleppointer to the structure to update
Returns
0 if no error detected, > 0 otherwise

◆ cs_param_saddle_set_solver()

int cs_param_saddle_set_solver ( const char *  keyval,
cs_param_saddle_t saddlep 
)

Set the type of solver to apply for this saddle-point system.

Parameters
[in]keyvalvalue of the key for the preconditioner
[in,out]saddleppointer to the structure to update
Returns
0 if no error detected, > 0 otherwise

◆ cs_param_saddle_set_solver_class()

int cs_param_saddle_set_solver_class ( const char *  keyval,
cs_param_saddle_t saddlep 
)

Set the class of solver to apply for this saddle-point system.

Parameters
[in]keyvalvalue of the key for the preconditioner
[in,out]saddleppointer to the structure to update
Returns
0 if no error detected, > 0 otherwise

◆ cs_param_saddle_try_init_schur_sles_param()

void cs_param_saddle_try_init_schur_sles_param ( cs_param_saddle_t saddlep)

Initialize a cs_param_sles_t structure for the Schur approximation nested inside a cs_param_saddle_t structure. By default, this member is not allocated. Do nothing if the related structure is already allocated.

Parameters
[in,out]saddleppointer to the structure to update