8.3
general documentation
cs_param_sles.h File Reference

Structure and routines handling the SLES ((Sparse Linear Equation Solver) settings stored inside a cs_param_sles_t structure. More...

#include "cs_param_amg.h"
#include "cs_param_mumps.h"
#include "cs_param_types.h"
+ Include dependency graph for cs_param_sles.h:

Go to the source code of this file.

Data Structures

struct  cs_param_sles_t
 Structure storing all metadata related to the resolution of a linear system with an iterative solver. More...
 

Functions

cs_param_sles_tcs_param_sles_create (int field_id, const char *system_name)
 Create a cs_param_sles_t structure and assign a default settings. More...
 
void cs_param_sles_free (cs_param_sles_t **p_slesp)
 Free a cs_param_sles_t structure. More...
 
void cs_param_sles_log (cs_param_sles_t *slesp)
 Log information related to the linear settings stored in the structure. More...
 
void cs_param_sles_copy_from (const cs_param_sles_t *src, cs_param_sles_t *dst)
 Copy a cs_param_sles_t structure from src to dst. More...
 
int cs_param_sles_set_solver (const char *keyval, cs_param_sles_t *slesp)
 Set the solver associated to this SLES from its keyval. More...
 
int cs_param_sles_set_precond (const char *keyval, cs_param_sles_t *slesp)
 Set the preconditioner associated to this SLES from its keyval. More...
 
int cs_param_sles_set_solver_class (const char *keyval, cs_param_sles_t *slesp)
 Set the class of solvers associated to this SLES from its keyval Common choices are "petsc", "hypre" or "mumps". More...
 
int cs_param_sles_set_amg_type (const char *keyval, cs_param_sles_t *slesp)
 Set the type of algebraic multigrid (AMG) associated to this SLES from its keyval. More...
 
int cs_param_sles_set_precond_block_type (const char *keyval, cs_param_sles_t *slesp)
 Set the type of block preconditioner associated to this SLES from its keyval. More...
 
void cs_param_sles_set_cvg_param (cs_param_sles_t *slesp, double rtol, double atol, double dtol, int max_iter)
 Set the convergence criteria for the given SLES parameters. One can use the parameter value CS_CDO_KEEP_DEFAULT to let unchange the settings. More...
 
void cs_param_sles_amg_inhouse_reset (cs_param_sles_t *slesp, bool used_as_solver, bool used_as_k_cycle)
 Allocate and initialize a new context structure for the in-house AMG settings. More...
 
void cs_param_sles_amg_inhouse (cs_param_sles_t *slesp, int n_down_iter, cs_param_amg_inhouse_solver_t down_smoother, int down_poly_deg, int n_up_iter, cs_param_amg_inhouse_solver_t up_smoother, int up_poly_deg, cs_param_amg_inhouse_solver_t coarse_solver, int coarse_poly_deg, cs_param_amg_inhouse_coarsen_t coarsen_algo, int aggreg_limit)
 Set the main members of a cs_param_amg_inhouse_t structure. This structure is not reset before applying the given parameter. Please call cs_param_sles_amg_inhouse_reset to do so. More...
 
void cs_param_sles_amg_inhouse_advanced (cs_param_sles_t *slesp, int max_levels, cs_gnum_t min_n_g_rows, double p0p1_relax, int coarse_max_iter, double coarse_rtol_mult)
 Set the members of a cs_param_amg_inhouse_t structure used in advanced settings. CS_CDO_KEEP_DEFAULT can be used to let unchange the parameter value. More...
 
void cs_param_sles_boomeramg_reset (cs_param_sles_t *slesp)
 Allocate and initialize a new context structure for the boomerAMG settings. More...
 
void cs_param_sles_boomeramg (cs_param_sles_t *slesp, int n_down_iter, cs_param_amg_boomer_smoother_t down_smoother, int n_up_iter, cs_param_amg_boomer_smoother_t up_smoother, cs_param_amg_boomer_smoother_t coarse_solver, cs_param_amg_boomer_coarsen_algo_t coarsen_algo)
 Set the main members of a cs_param_amg_boomer_t structure. This structure is allocated and initialized and then one sets the main given parameters. Please refer to the HYPRE user guide for more details about the following options. More...
 
void cs_param_sles_boomeramg_advanced (cs_param_sles_t *slesp, double strong_thr, cs_param_amg_boomer_interp_algo_t interp_algo, int p_max, int n_agg_lv, int n_agg_paths)
 Set the members of a cs_param_amg_boomer_t structure used in advanced settings. This structure is allocated if needed. Other members are kept to their values. Please refer to the HYPRE user guide for more details about the following options. More...
 
void cs_param_sles_gamg_reset (cs_param_sles_t *slesp)
 Allocate and initialize a new context structure for the GAMG settings in PETSc. More...
 
void cs_param_sles_gamg (cs_param_sles_t *slesp, int n_down_iter, cs_param_amg_gamg_smoother_t down_smoother, int n_up_iter, cs_param_amg_gamg_smoother_t up_smoother, cs_param_amg_gamg_coarse_solver_t coarse_solver)
 Set the main members of a cs_param_amg_gamg_t structure. This structure is allocated, initialized by default and then one sets the main given parameters. Please refer to the PETSc user guide for more details about the following options. More...
 
void cs_param_sles_gamg_advanced (cs_param_sles_t *slesp, double threshold, int n_agg_lv, bool use_sq_grph, int n_smooth_agg)
 Set the members of a cs_param_amg_gamg_t structure used in advanced settings. This structure is allocated if needed. Other members are kept to their values. Please refer to the PETSc user guide for more details about the following options. More...
 
void cs_param_sles_hmg_reset (cs_param_sles_t *slesp)
 Allocate and initialize a new context structure for the HMG settings in PETSc. More...
 
void cs_param_sles_hmg (cs_param_sles_t *slesp, int n_down_iter, cs_param_amg_gamg_smoother_t down_smoother, int n_up_iter, cs_param_amg_gamg_smoother_t up_smoother, cs_param_amg_gamg_coarse_solver_t coarse_solver)
 Set the main members of a cs_param_amg_hmg_t structure. This structure is allocated, initialized by default and then one sets the main given parameters. Please refer to the PETSc user guide for more details about the following options. More...
 
void cs_param_sles_hmg_advanced (cs_param_sles_t *slesp, bool use_boomer_coarsening, bool reuse_interpolation, bool subspace_coarsening)
 Set the members of a cs_param_amg_gamg_t structure used in advanced settings. This structure is allocated if needed. Other members are kept to their values. Please refer to the PETSc user guide for more details about the following options. More...
 
void cs_param_sles_mumps_reset (cs_param_sles_t *slesp)
 Allocate and initialize a new context structure for the MUMPS settings. More...
 
void cs_param_sles_mumps (cs_param_sles_t *slesp, bool is_single, cs_param_mumps_facto_type_t facto_type)
 Set the main members of a cs_param_mumps_t structure. This structure is allocated and initialized with default settings if needed. If the structure exists already, then advanced members are kept to their values. More...
 
void cs_param_sles_mumps_advanced (cs_param_sles_t *slesp, cs_param_mumps_analysis_algo_t analysis_algo, int block_analysis, bool keep_ordering, double mem_coef, double blr_threshold, int ir_steps, cs_param_mumps_memory_usage_t mem_usage, bool advanced_optim)
 Set the members related to an advanced settings of a cs_param_mumps_t structure. This structure is allocated and initialized if needed. Please refer to the MUMPS user guide for more details about the following advanced options. More...
 
bool cs_param_sles_hypre_from_petsc (void)
 Check the availability of Hypre solvers from the PETSc library. More...
 
cs_param_solver_class_t cs_param_sles_check_class (cs_param_solver_class_t wanted_class)
 Check the availability of a solver library and return the requested one if this is possible or an alternative or CS_PARAM_N_SOLVER_CLASSES if no alternative is available. More...
 

Detailed Description

Structure and routines handling the SLES ((Sparse Linear Equation Solver) settings stored inside a cs_param_sles_t structure.

Function Documentation

◆ cs_param_sles_amg_inhouse()

void cs_param_sles_amg_inhouse ( cs_param_sles_t slesp,
int  n_down_iter,
cs_param_amg_inhouse_solver_t  down_smoother,
int  down_poly_deg,
int  n_up_iter,
cs_param_amg_inhouse_solver_t  up_smoother,
int  up_poly_deg,
cs_param_amg_inhouse_solver_t  coarse_solver,
int  coarse_poly_deg,
cs_param_amg_inhouse_coarsen_t  coarsen_algo,
int  aggreg_limit 
)

Set the main members of a cs_param_amg_inhouse_t structure. This structure is not reset before applying the given parameter. Please call cs_param_sles_amg_inhouse_reset to do so.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]n_down_iternumber of smoothing steps for the down cycle
[in]down_smoothertype of smoother for the down cycle
[in]down_poly_degpoly. degree for the down smoother
[in]n_up_iternumber of smoothing steps for the up cycle
[in]up_smoothertype of smoother for the up cycle
[in]up_poly_degpoly. degree for the up smoother
[in]coarse_solversolver at the coarsest level
[in]coarse_poly_degpoly. degree for the coarse solver
[in]coarsen_algotype of algorithm to coarsen grids
[in]aggreg_limitlimit for the aggregation process

◆ cs_param_sles_amg_inhouse_advanced()

void cs_param_sles_amg_inhouse_advanced ( cs_param_sles_t slesp,
int  max_levels,
cs_gnum_t  min_n_g_rows,
double  p0p1_relax,
int  coarse_max_iter,
double  coarse_rtol_mult 
)

Set the members of a cs_param_amg_inhouse_t structure used in advanced settings. CS_CDO_KEEP_DEFAULT can be used to let unchange the parameter value.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]max_levelsmax. number of levels
[in]min_n_g_rowsdo not coarsen anymore below this number
[in]p0p1_relaxp0/p1 relaxation parameter
[in]coarse_max_itermax. number of iter. for the coarse solver
[in]coarse_rtol_multmax. number of iter. for the coarse solver

◆ cs_param_sles_amg_inhouse_reset()

void cs_param_sles_amg_inhouse_reset ( cs_param_sles_t slesp,
bool  used_as_solver,
bool  used_as_k_cycle 
)

Allocate and initialize a new context structure for the in-house AMG settings.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]used_as_solvertrue or false
[in]used_as_k_cycletrue or false

◆ cs_param_sles_boomeramg()

void cs_param_sles_boomeramg ( cs_param_sles_t slesp,
int  n_down_iter,
cs_param_amg_boomer_smoother_t  down_smoother,
int  n_up_iter,
cs_param_amg_boomer_smoother_t  up_smoother,
cs_param_amg_boomer_smoother_t  coarse_solver,
cs_param_amg_boomer_coarsen_algo_t  coarsen_algo 
)

Set the main members of a cs_param_amg_boomer_t structure. This structure is allocated and initialized and then one sets the main given parameters. Please refer to the HYPRE user guide for more details about the following options.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]n_down_iternumber of smoothing steps for the down cycle
[in]down_smoothertype of smoother for the down cycle
[in]n_up_iternumber of smoothing steps for the up cycle
[in]up_smoothertype of smoother for the up cycle
[in]coarse_solversolver at the coarsest level
[in]coarsen_algotype of algorithm to coarsen grids

◆ cs_param_sles_boomeramg_advanced()

void cs_param_sles_boomeramg_advanced ( cs_param_sles_t slesp,
double  strong_thr,
cs_param_amg_boomer_interp_algo_t  interp_algo,
int  p_max,
int  n_agg_lv,
int  n_agg_paths 
)

Set the members of a cs_param_amg_boomer_t structure used in advanced settings. This structure is allocated if needed. Other members are kept to their values. Please refer to the HYPRE user guide for more details about the following options.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]strong_thrvalue of the strong threshold (coarsening)
[in]interp_algoalgorithm used for the interpolation
[in]p_maxmax number of elements per row (interp)
[in]n_agg_lvaggressive coarsening (number of levels)
[in]n_agg_pathsaggressive coarsening (number of paths)

◆ cs_param_sles_boomeramg_reset()

void cs_param_sles_boomeramg_reset ( cs_param_sles_t slesp)

Allocate and initialize a new context structure for the boomerAMG settings.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure

◆ cs_param_sles_check_class()

cs_param_solver_class_t cs_param_sles_check_class ( cs_param_solver_class_t  wanted_class)

Check the availability of a solver library and return the requested one if this is possible or an alternative or CS_PARAM_N_SOLVER_CLASSES if no alternative is available.

Parameters
[in]wanted_classrequested class of solvers
Returns
the available solver class related to the requested class

◆ cs_param_sles_copy_from()

void cs_param_sles_copy_from ( const cs_param_sles_t src,
cs_param_sles_t dst 
)

Copy a cs_param_sles_t structure from src to dst.

Parameters
[in]srcreference cs_param_sles_t structure to copy
[in,out]dstcopy of the reference at exit

◆ cs_param_sles_create()

cs_param_sles_t * cs_param_sles_create ( int  field_id,
const char *  system_name 
)

Create a cs_param_sles_t structure and assign a default settings.

Parameters
[in]field_idid related to to the variable field or -1
[in]system_namename of the system to solve or NULL
Returns
a pointer to a cs_param_sles_t stucture
Parameters
[in]field_idid related to to the variable field or -1
[in]system_namename of the system to solve or nullptr
Returns
a pointer to a cs_param_sles_t stucture

◆ cs_param_sles_free()

void cs_param_sles_free ( cs_param_sles_t **  p_slesp)

Free a cs_param_sles_t structure.

Parameters
[in,out]p_slesppointer to a cs_param_sles_t structure to free

◆ cs_param_sles_gamg()

void cs_param_sles_gamg ( cs_param_sles_t slesp,
int  n_down_iter,
cs_param_amg_gamg_smoother_t  down_smoother,
int  n_up_iter,
cs_param_amg_gamg_smoother_t  up_smoother,
cs_param_amg_gamg_coarse_solver_t  coarse_solver 
)

Set the main members of a cs_param_amg_gamg_t structure. This structure is allocated, initialized by default and then one sets the main given parameters. Please refer to the PETSc user guide for more details about the following options.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]n_down_iternumber of smoothing steps for the down cycle
[in]down_smoothertype of smoother for the down cycle
[in]n_up_iternumber of smoothing steps for the up cycle
[in]up_smoothertype of smoother for the up cycle
[in]coarse_solversolver at the coarsest level

◆ cs_param_sles_gamg_advanced()

void cs_param_sles_gamg_advanced ( cs_param_sles_t slesp,
double  threshold,
int  n_agg_lv,
bool  use_sq_grph,
int  n_smooth_agg 
)

Set the members of a cs_param_amg_gamg_t structure used in advanced settings. This structure is allocated if needed. Other members are kept to their values. Please refer to the PETSc user guide for more details about the following options.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]thresholdvalue of the coarsening threshold
[in]n_agg_lvaggressive coarsening (number of levels)
[in]use_sq_grphuse previous square graph for aggressive coa.
[in]n_smooth_aggsmooth aggregation (number of sweeps)

◆ cs_param_sles_gamg_reset()

void cs_param_sles_gamg_reset ( cs_param_sles_t slesp)

Allocate and initialize a new context structure for the GAMG settings in PETSc.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure

◆ cs_param_sles_hmg()

void cs_param_sles_hmg ( cs_param_sles_t slesp,
int  n_down_iter,
cs_param_amg_gamg_smoother_t  down_smoother,
int  n_up_iter,
cs_param_amg_gamg_smoother_t  up_smoother,
cs_param_amg_gamg_coarse_solver_t  coarse_solver 
)

Set the main members of a cs_param_amg_hmg_t structure. This structure is allocated, initialized by default and then one sets the main given parameters. Please refer to the PETSc user guide for more details about the following options.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]n_down_iternumber of smoothing steps for the down cycle
[in]down_smoothertype of smoother for the down cycle
[in]n_up_iternumber of smoothing steps for the up cycle
[in]up_smoothertype of smoother for the up cycle
[in]coarse_solversolver at the coarsest level

◆ cs_param_sles_hmg_advanced()

void cs_param_sles_hmg_advanced ( cs_param_sles_t slesp,
bool  use_boomer_coarsening,
bool  reuse_interpolation,
bool  subspace_coarsening 
)

Set the members of a cs_param_amg_gamg_t structure used in advanced settings. This structure is allocated if needed. Other members are kept to their values. Please refer to the PETSc user guide for more details about the following options.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]use_boomer_coarseningcoarsening done by HYPRE ?
[in]reuse_interpolationreuse interpolation for new mat. values
[in]subspace_coarseninguse coarsening of submatrices

◆ cs_param_sles_hmg_reset()

void cs_param_sles_hmg_reset ( cs_param_sles_t slesp)

Allocate and initialize a new context structure for the HMG settings in PETSc.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure

◆ cs_param_sles_hypre_from_petsc()

bool cs_param_sles_hypre_from_petsc ( void  )

Check the availability of Hypre solvers from the PETSc library.

Returns
return true or false

◆ cs_param_sles_log()

void cs_param_sles_log ( cs_param_sles_t slesp)

Log information related to the linear settings stored in the structure.

Parameters
[in]slesppointer to a cs_param_sles_log

◆ cs_param_sles_mumps()

void cs_param_sles_mumps ( cs_param_sles_t slesp,
bool  is_single,
cs_param_mumps_facto_type_t  facto_type 
)

Set the main members of a cs_param_mumps_t structure. This structure is allocated and initialized with default settings if needed. If the structure exists already, then advanced members are kept to their values.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]is_singlesingle-precision or double-precision
[in]facto_typetype of factorization to consider

◆ cs_param_sles_mumps_advanced()

void cs_param_sles_mumps_advanced ( cs_param_sles_t slesp,
cs_param_mumps_analysis_algo_t  analysis_algo,
int  block_analysis,
bool  keep_ordering,
double  mem_coef,
double  blr_threshold,
int  ir_steps,
cs_param_mumps_memory_usage_t  mem_usage,
bool  advanced_optim 
)

Set the members related to an advanced settings of a cs_param_mumps_t structure. This structure is allocated and initialized if needed. Please refer to the MUMPS user guide for more details about the following advanced options.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]analysis_algoalgorithm used for the analysis step
[in]block_analysis> 1: fixed block size; otherwise do nothing
[in]keep_orderingtrue: keep the initial ordering to save time
[in]mem_coefpercentage increase in the memory workspace
[in]blr_thresholdAccuracy in BLR compression (0: not used)
[in]ir_steps0: No, otherwise the number of iterations
[in]mem_usagestrategy to adopt for the memory usage
[in]advanced_optimactivate advanced optimization (MPI/openMP)

◆ cs_param_sles_mumps_reset()

void cs_param_sles_mumps_reset ( cs_param_sles_t slesp)

Allocate and initialize a new context structure for the MUMPS settings.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure

◆ cs_param_sles_set_amg_type()

int cs_param_sles_set_amg_type ( const char *  keyval,
cs_param_sles_t slesp 
)

Set the type of algebraic multigrid (AMG) associated to this SLES from its keyval.

Parameters
[in]keyvalvalue of the key
[in,out]slesppointer to a cs_param_sles_t structure
Returns
an error code (> 0) or 0 if there is no error

◆ cs_param_sles_set_cvg_param()

void cs_param_sles_set_cvg_param ( cs_param_sles_t slesp,
double  rtol,
double  atol,
double  dtol,
int  max_iter 
)

Set the convergence criteria for the given SLES parameters. One can use the parameter value CS_CDO_KEEP_DEFAULT to let unchange the settings.

Parameters
[in,out]slesppointer to a cs_param_sles_t structure
[in]rtolrelative tolerance
[in]atolabsolute tolerance
[in]dtoldivergence tolerance
[in]max_itermax. number of iterations

◆ cs_param_sles_set_precond()

int cs_param_sles_set_precond ( const char *  keyval,
cs_param_sles_t slesp 
)

Set the preconditioner associated to this SLES from its keyval.

Parameters
[in]keyvalvalue of the key
[in,out]slesppointer to a cs_param_sles_t structure
Returns
an error code (> 0) or 0 if there is no error

◆ cs_param_sles_set_precond_block_type()

int cs_param_sles_set_precond_block_type ( const char *  keyval,
cs_param_sles_t slesp 
)

Set the type of block preconditioner associated to this SLES from its keyval.

Parameters
[in]keyvalvalue of the key
[in,out]slesppointer to a cs_param_sles_t structure
Returns
an error code (> 0) or 0 if there is no error

◆ cs_param_sles_set_solver()

int cs_param_sles_set_solver ( const char *  keyval,
cs_param_sles_t slesp 
)

Set the solver associated to this SLES from its keyval.

Parameters
[in]keyvalvalue of the key
[in,out]slesppointer to a cs_param_sles_t structure
Returns
an error code (> 0) or 0 if there is no error

◆ cs_param_sles_set_solver_class()

int cs_param_sles_set_solver_class ( const char *  keyval,
cs_param_sles_t slesp 
)

Set the class of solvers associated to this SLES from its keyval Common choices are "petsc", "hypre" or "mumps".

Parameters
[in]keyvalvalue of the key
[in,out]slesppointer to a cs_param_sles_t structure
Returns
an error code (> 0) or 0 if there is no error