7.3
general documentation
cs_sles_it.h File Reference
#include "cs_base.h"
#include "cs_matrix.h"
#include "cs_sles.h"
#include "cs_sles_pc.h"
+ Include dependency graph for cs_sles_it.h:

Go to the source code of this file.

Typedefs

typedef struct _cs_sles_it_t cs_sles_it_t
 
typedef struct _cs_sles_it_convergence_t cs_sles_it_convergence_t
 

Enumerations

enum  cs_sles_it_type_t {
  CS_SLES_PCG, CS_SLES_FCG, CS_SLES_IPCG, CS_SLES_JACOBI,
  CS_SLES_BICGSTAB, CS_SLES_BICGSTAB2, CS_SLES_GCR, CS_SLES_GMRES,
  CS_SLES_P_GAUSS_SEIDEL, CS_SLES_P_SYM_GAUSS_SEIDEL, CS_SLES_PCR3, CS_SLES_USER_DEFINED,
  CS_SLES_N_IT_TYPES, CS_SLES_TS_F_GAUSS_SEIDEL, CS_SLES_TS_B_GAUSS_SEIDEL, CS_SLES_N_SMOOTHER_TYPES
}
 

Functions

cs_sles_convergence_state_t cs_user_sles_it_solver (cs_sles_it_t *c, const cs_matrix_t *a, cs_lnum_t diag_block_size, cs_sles_it_convergence_t *convergence, const cs_real_t *rhs, cs_real_t *restrict vx, size_t aux_size, void *aux_vectors)
 
cs_sles_it_tcs_sles_it_define (int f_id, const char *name, cs_sles_it_type_t solver_type, int poly_degree, int n_max_iter)
 Define and associate an iterative sparse linear system solver for a given field or equation name. More...
 
cs_sles_it_tcs_sles_it_create (cs_sles_it_type_t solver_type, int poly_degree, int n_max_iter, bool update_stats)
 Create iterative sparse linear system solver info and context. More...
 
void cs_sles_it_destroy (void **context)
 Destroy iterative sparse linear system solver info and context. More...
 
void * cs_sles_it_copy (const void *context)
 Create iterative sparse linear system solver info and context based on existing info and context. More...
 
void cs_sles_it_setup (void *context, const char *name, const cs_matrix_t *a, int verbosity)
 Setup iterative sparse linear equation solver. More...
 
cs_sles_convergence_state_t cs_sles_it_solve (void *context, const char *name, const cs_matrix_t *a, int verbosity, double precision, double r_norm, int *n_iter, double *residue, const cs_real_t *rhs, cs_real_t *vx, size_t aux_size, void *aux_vectors)
 Call iterative sparse linear equation solver. More...
 
void cs_sles_it_free (void *context)
 Free iterative sparse linear equation solver setup context. More...
 
void cs_sles_it_log (const void *context, cs_log_t log_type)
 Log sparse linear equation solver info. More...
 
cs_sles_it_type_t cs_sles_it_get_type (const cs_sles_it_t *context)
 Return iterative solver type. More...
 
double cs_sles_it_get_last_initial_residue (const cs_sles_it_t *context)
 Return the initial residue for the previous solve operation with a solver. More...
 
cs_sles_pc_tcs_sles_it_get_pc (cs_sles_it_t *context)
 Return a preconditioner context for an iterative sparse linear equation solver. More...
 
void cs_sles_it_transfer_pc (cs_sles_it_t *context, cs_sles_pc_t **pc)
 Assign a preconditioner to an iterative sparse linear equation solver, transfering its ownership to to solver context. More...
 
void cs_sles_it_transfer_parameters (const cs_sles_it_t *src, cs_sles_it_t *dest)
 Copy options from one iterative sparse linear system solver info and context to another. More...
 
void cs_sles_it_set_shareable (cs_sles_it_t *context, const cs_sles_it_t *shareable)
 Associate a similar info and context object with which some setup data may be shared. More...
 
void cs_sles_it_assign_order (cs_sles_it_t *context, cs_lnum_t **order)
 Assign ordering to iterative solver. More...
 
double cs_sles_it_get_breakdown_threshold (void)
 Retrieve the threshold value under which a breakdown happens in solvers like BiCGStab or BiCGStab2. More...
 
void cs_sles_it_set_breakdown_threshold (double threshold)
 Define the threshold value under which a breakdown happens in solvers like BiCGStab or BiCGStab2. More...
 
void cs_sles_it_set_fallback_threshold (cs_sles_it_t *context, cs_sles_convergence_state_t threshold)
 Define convergence level under which the fallback to another solver may be used if applicable. More...
 
void cs_sles_it_set_restart_interval (cs_sles_it_t *context, int interval)
 Define the number of iterations to be done before restarting the solver. Useful only for GCR or GMRES algorithms. More...
 
void cs_sles_it_set_n_max_iter (cs_sles_it_t *context, int n_max_iter)
 Define the max. number of iterations before stopping the algorithm. More...
 
cs_lnum_t cs_sles_it_get_pcg_single_reduction (void)
 Query mean number of rows under which Conjugate Gradient algorithm uses the single-reduction variant. More...
 
void cs_sles_it_set_pcg_single_reduction (cs_lnum_t threshold)
 Set mean number of rows under which Conjugate Gradient algorithm should use the single-reduction variant. More...
 
void cs_sles_it_log_parallel_options (void)
 Log the current global settings relative to parallelism. More...
 
bool cs_sles_it_error_post_and_abort (cs_sles_t *sles, cs_sles_convergence_state_t state, const cs_matrix_t *a, const cs_real_t *rhs, cs_real_t *vx)
 Error handler for iterative sparse linear equation solver. More...
 
void cs_sles_it_set_plot_options (cs_sles_it_t *context, const char *base_name, bool use_iteration)
 Set plotting options for an iterative sparse linear equation solver. More...
 
cs_sles_convergence_state_t cs_sles_it_convergence_test (cs_sles_it_t *c, unsigned n_iter, double residue, cs_sles_it_convergence_t *convergence)
 

Variables

const char * cs_sles_it_type_name []
 

Typedef Documentation

◆ cs_sles_it_convergence_t

typedef struct _cs_sles_it_convergence_t cs_sles_it_convergence_t

◆ cs_sles_it_t

typedef struct _cs_sles_it_t cs_sles_it_t

Enumeration Type Documentation

◆ cs_sles_it_type_t

Enumerator
CS_SLES_PCG 

Preconditioned conjugate gradient

CS_SLES_FCG 

Preconditions flexible conjugate gradient, described in [10]

CS_SLES_IPCG 

Preconditions inexact conjugate gradient

CS_SLES_JACOBI 

Jacobi

CS_SLES_BICGSTAB 

Preconditioned BiCGstab (biconjugate gradient stabilized)

CS_SLES_BICGSTAB2 

Preconditioned BiCGstab2

CS_SLES_GCR 

Generalized conjugate residual

CS_SLES_GMRES 

Preconditioned GMRES (generalized minimal residual)

CS_SLES_P_GAUSS_SEIDEL 

Process-local Gauss-Seidel

CS_SLES_P_SYM_GAUSS_SEIDEL 

Process-local symmetric Gauss-Seidel

CS_SLES_PCR3 

3-layer conjugate residual

CS_SLES_USER_DEFINED 

User-defined iterative solver

CS_SLES_N_IT_TYPES 

Number of resolution algorithms excluding smoother only

CS_SLES_TS_F_GAUSS_SEIDEL 

Truncated forward Gauss-Seidel smoother

CS_SLES_TS_B_GAUSS_SEIDEL 

Truncated backward Gauss-Seidel smoother

CS_SLES_N_SMOOTHER_TYPES 

Number of resolution algorithms including smoother only

Function Documentation

◆ cs_sles_it_assign_order()

void cs_sles_it_assign_order ( cs_sles_it_t context,
cs_lnum_t **  order 
)

Assign ordering to iterative solver.

The solver context takes ownership of the order array (i.e. it will handle its later deallocation).

This is useful only for Process-local Gauss-Seidel.

Parameters
[in,out]contextpointer to iterative solver info and context
[in,out]orderpointer to ordering array

◆ cs_sles_it_convergence_test()

cs_sles_convergence_state_t cs_sles_it_convergence_test ( cs_sles_it_t c,
unsigned  n_iter,
double  residue,
cs_sles_it_convergence_t convergence 
)

◆ cs_sles_it_copy()

void* cs_sles_it_copy ( const void *  context)

Create iterative sparse linear system solver info and context based on existing info and context.

Parameters
[in]contextpointer to reference info and context (actual type: cs_sles_it_t *)
Returns
pointer to newly created solver info object. (actual type: cs_sles_it_t *)

◆ cs_sles_it_create()

cs_sles_it_t* cs_sles_it_create ( cs_sles_it_type_t  solver_type,
int  poly_degree,
int  n_max_iter,
bool  update_stats 
)

Create iterative sparse linear system solver info and context.

Parameters
[in]solver_typetype of solver (PCG, Jacobi, ...)
[in]poly_degreepreconditioning polynomial degree (0: diagonal; -1: non-preconditioned; see Iterative linear solvers. for details)
[in]n_max_itermaximum number of iterations
[in]update_statsautomatic solver statistics indicator
Returns
pointer to newly created solver info object.

◆ cs_sles_it_define()

cs_sles_it_t* cs_sles_it_define ( int  f_id,
const char *  name,
cs_sles_it_type_t  solver_type,
int  poly_degree,
int  n_max_iter 
)

Define and associate an iterative sparse linear system solver for a given field or equation name.

If this system did not previously exist, it is added to the list of "known" systems. Otherwise, its definition is replaced by the one defined here.

This is a utility function: if finer control is needed, see cs_sles_define and cs_sles_it_create.

Note that this function returns a pointer directly to the iterative solver management structure. This may be used to set further options, for example using cs_sles_it_set_plot_options. If needed, cs_sles_find may be used to obtain a pointer to the matching cs_sles_t container.

Parameters
[in]f_idassociated field id, or < 0
[in]nameassociated name if f_id < 0, or NULL
[in]solver_typetype of solver (PCG, Jacobi, ...)
[in]poly_degreepreconditioning polynomial degree (0: diagonal; -1: non-preconditioned)
[in]n_max_itermaximum number of iterations
Returns
pointer to newly created iterative solver info object.

◆ cs_sles_it_destroy()

void cs_sles_it_destroy ( void **  context)

Destroy iterative sparse linear system solver info and context.

Parameters
[in,out]contextpointer to iterative solver info and context (actual type: cs_sles_it_t **)

◆ cs_sles_it_error_post_and_abort()

bool cs_sles_it_error_post_and_abort ( cs_sles_t sles,
cs_sles_convergence_state_t  state,
const cs_matrix_t a,
const cs_real_t rhs,
cs_real_t vx 
)

Error handler for iterative sparse linear equation solver.

In case of divergence or breakdown, this error handler outputs postprocessing data to assist debugging, then aborts the run. It does nothing in case the maximum iteration count is reached.

Parameters
[in,out]slespointer to solver object
[in]stateconvergence state
[in]amatrix
[in]rhsright hand side
[in,out]vxsystem solution
Returns
false (do not attempt new solve)

◆ cs_sles_it_free()

void cs_sles_it_free ( void *  context)

Free iterative sparse linear equation solver setup context.

This function frees resolution-related data, such as buffers and preconditioning but does not free the whole context, as info used for logging (especially performance data) is maintained.

Parameters
[in,out]contextpointer to iterative solver info and context (actual type: cs_sles_it_t *)

◆ cs_sles_it_get_breakdown_threshold()

double cs_sles_it_get_breakdown_threshold ( void  )

Retrieve the threshold value under which a breakdown happens in solvers like BiCGStab or BiCGStab2.

Returns
the value of the threshold

◆ cs_sles_it_get_last_initial_residue()

double cs_sles_it_get_last_initial_residue ( const cs_sles_it_t context)

Return the initial residue for the previous solve operation with a solver.

This is useful for convergence tests when this solver is used as a preconditioning smoother.

This operation is only valid between calls to cs_sles_it_setup (or cs_sles_it_solve) and cs_sles_it_free. It returns -1 otherwise.

Parameters
[in]contextpointer to iterative solver info and context
Returns
initial residue from last call to cs_sles_solve with this solver

◆ cs_sles_it_get_pc()

cs_sles_pc_t* cs_sles_it_get_pc ( cs_sles_it_t context)

Return a preconditioner context for an iterative sparse linear equation solver.

This allows modifying parameters of a non default (Jacobi or polynomial) preconditioner.

Parameters
[in]contextpointer to iterative solver info and context
Returns
pointer to preconditoner context

◆ cs_sles_it_get_pcg_single_reduction()

cs_lnum_t cs_sles_it_get_pcg_single_reduction ( void  )

Query mean number of rows under which Conjugate Gradient algorithm uses the single-reduction variant.

The single-reduction variant requires only one parallel sum per iteration (instead of 2), at the cost of additional vector operations, so it tends to be more expensive when the number of matrix rows per MPI rank is high, then becomes cheaper when the MPI latency cost becomes more significant.

This option is ignored for non-parallel runs, so 0 is returned.

Returns
mean number of rows per active rank under which the single-reduction variant will be used

◆ cs_sles_it_get_type()

cs_sles_it_type_t cs_sles_it_get_type ( const cs_sles_it_t context)

Return iterative solver type.

Parameters
[in]contextpointer to iterative solver info and context
Returns
selected solver type

◆ cs_sles_it_log()

void cs_sles_it_log ( const void *  context,
cs_log_t  log_type 
)

Log sparse linear equation solver info.

Parameters
[in]contextpointer to iterative solver info and context (actual type: cs_sles_it_t *)
[in]log_typelog type

◆ cs_sles_it_log_parallel_options()

void cs_sles_it_log_parallel_options ( void  )

Log the current global settings relative to parallelism.

◆ cs_sles_it_set_breakdown_threshold()

void cs_sles_it_set_breakdown_threshold ( double  threshold)

Define the threshold value under which a breakdown happens in solvers like BiCGStab or BiCGStab2.

Parameters
[in]thresholdvalue of the threshold

◆ cs_sles_it_set_fallback_threshold()

void cs_sles_it_set_fallback_threshold ( cs_sles_it_t context,
cs_sles_convergence_state_t  threshold 
)

Define convergence level under which the fallback to another solver may be used if applicable.

Currently, this mechanism is only by default used for BiCGstab and 3-layer conjugate residual solvers with scalar matrices, which may fall back to a preconditioned GMRES solver. For those solvers, the default threshold is CS_SLES_BREAKDOWN, meaning that divergence (but not breakdown) will lead to the use of the fallback mechanism.

Parameters
[in,out]contextpointer to iterative solver info and context
[in]thresholdconvergence level under which fallback is used

◆ cs_sles_it_set_n_max_iter()

void cs_sles_it_set_n_max_iter ( cs_sles_it_t context,
int  n_max_iter 
)

Define the max. number of iterations before stopping the algorithm.

Parameters
[in,out]contextpointer to iterative solver info and context
[in]n_max_itermax. number of iterations

◆ cs_sles_it_set_pcg_single_reduction()

void cs_sles_it_set_pcg_single_reduction ( cs_lnum_t  threshold)

Set mean number of rows under which Conjugate Gradient algorithm should use the single-reduction variant.

The single-reduction variant requires only one parallel sum per iteration (instead of 2), at the cost of additional vector operations, so it tends to be more expensive when the number of matrix rows per MPI rank is high, then becomes cheaper when the MPI latency cost becomes more significant.

This option is ignored for non-parallel runs.

Parameters
[in]thresholdmean number of rows per active rank under which the single-reduction variant will be used

◆ cs_sles_it_set_plot_options()

void cs_sles_it_set_plot_options ( cs_sles_it_t context,
const char *  base_name,
bool  use_iteration 
)

Set plotting options for an iterative sparse linear equation solver.

Parameters
[in,out]contextpointer to iterative solver info and context
[in]base_namebase plot name to activate, NULL otherwise
[in]use_iterationif true, use iteration as time stamp otherwise, use wall clock time

◆ cs_sles_it_set_restart_interval()

void cs_sles_it_set_restart_interval ( cs_sles_it_t context,
int  interval 
)

Define the number of iterations to be done before restarting the solver. Useful only for GCR or GMRES algorithms.

Parameters
[in,out]contextpointer to iterative solver info and context
[in]intervalconvergence level under which fallback is used

◆ cs_sles_it_set_shareable()

void cs_sles_it_set_shareable ( cs_sles_it_t context,
const cs_sles_it_t shareable 
)

Associate a similar info and context object with which some setup data may be shared.

This is especially useful for sharing preconditioning data between similar solver contexts (for example ascending and descending multigrid smoothers based on the same matrix).

For preconditioning data to be effectively shared, cs_sles_it_setup (or cs_sles_it_solve) must be called on shareable before being called on context (without cs_sles_it_free being called in between, of course).

It is the caller's responsibility to ensure the context is not used for a cs_sles_it_setup or cs_sles_it_solve operation after the shareable object has been destroyed (normally by cs_sles_it_destroy).

Parameters
[in,out]contextpointer to iterative solver info and context
[in]shareablepointer to iterative solver info and context whose context may be shared

◆ cs_sles_it_setup()

void cs_sles_it_setup ( void *  context,
const char *  name,
const cs_matrix_t a,
int  verbosity 
)

Setup iterative sparse linear equation solver.

Parameters
[in,out]contextpointer to iterative solver info and context (actual type: cs_sles_it_t *)
[in]namepointer to system name
[in]aassociated matrix
[in]verbosityassociated verbosity

◆ cs_sles_it_solve()

cs_sles_convergence_state_t cs_sles_it_solve ( void *  context,
const char *  name,
const cs_matrix_t a,
int  verbosity,
double  precision,
double  r_norm,
int *  n_iter,
double *  residue,
const cs_real_t rhs,
cs_real_t vx,
size_t  aux_size,
void *  aux_vectors 
)

Call iterative sparse linear equation solver.

Parameters
[in,out]contextpointer to iterative solver info and context (actual type: cs_sles_it_t *)
[in]namepointer to system name
[in]amatrix
[in]verbosityassociated verbosity
[in]precisionsolver precision
[in]r_normresidue normalization
[out]n_iternumber of "equivalent" iterations
[out]residueresidue
[in]rhsright hand side
[in,out]vxsystem solution
[in]aux_sizenumber of elements in aux_vectors (in bytes)
aux_vectorsoptional working area (internal allocation if NULL)
Returns
convergence state

◆ cs_sles_it_transfer_parameters()

void cs_sles_it_transfer_parameters ( const cs_sles_it_t src,
cs_sles_it_t dest 
)

Copy options from one iterative sparse linear system solver info and context to another.

Optional plotting contexts are shared between the source and destination contexts.

Preconditioner settings are to be handled separately.

Parameters
[in]srcpointer to source info and context
[in,out]destpointer to destination info and context

◆ cs_sles_it_transfer_pc()

void cs_sles_it_transfer_pc ( cs_sles_it_t context,
cs_sles_pc_t **  pc 
)

Assign a preconditioner to an iterative sparse linear equation solver, transfering its ownership to to solver context.

This allows assigning a non default (Jacobi or polynomial) preconditioner.

The input pointer is set to NULL to make it clear the caller does not own the preconditioner anymore, though the context can be accessed using cs_sles_it_get_pc.

Parameters
[in,out]contextpointer to iterative solver info and context
[in,out]pcpointer to preconditoner context

◆ cs_user_sles_it_solver()

cs_sles_convergence_state_t cs_user_sles_it_solver ( cs_sles_it_t c,
const cs_matrix_t a,
cs_lnum_t  diag_block_size,
cs_sles_it_convergence_t convergence,
const cs_real_t rhs,
cs_real_t *restrict  vx,
size_t  aux_size,
void *  aux_vectors 
)

Variable Documentation

◆ cs_sles_it_type_name

const char* cs_sles_it_type_name[]