#include "cs_defs.h"
#include <stdarg.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <assert.h>
#include <math.h>
#include "bft_mem.h"
#include "bft_error.h"
#include "bft_printf.h"
#include "cs_base.h"
#include "cs_field.h"
#include "cs_field_default.h"
#include "cs_grid.h"
#include "cs_halo.h"
#include "cs_internal_coupling.h"
#include "cs_log.h"
#include "cs_mesh.h"
#include "cs_mesh_adjacencies.h"
#include "cs_mesh_quantities.h"
#include "cs_matrix.h"
#include "cs_matrix_default.h"
#include "cs_matrix_util.h"
#include "cs_multigrid.h"
#include "cs_parameters.h"
#include "cs_sles.h"
#include "cs_sles_it.h"
#include "cs_sles_pc.h"
#include "cs_timer.h"
#include "cs_sles_default.h"
Functions | |
void | cs_sles_default (int f_id, const char *name, const cs_matrix_t *a) |
Default definition of a sparse linear equation solver. More... | |
void | cs_sles_default_setup (void) |
Default setup for sparse linear equation solver API. More... | |
void | cs_sles_default_finalize (void) |
Default finalization for sparse linear equation solver API. More... | |
int | cs_sles_default_get_verbosity (int f_id, const char *name) |
Return default verbosity associated to a field id, name couple. More... | |
cs_matrix_t * | cs_sles_default_get_matrix (int f_id, const char *name, cs_lnum_t db_size, cs_lnum_t eb_size, bool symmetric) |
Return pointer to matrix structure matching equation solve. More... | |
void | cs_sles_default_release_matrix (cs_matrix_t **m) |
void | cs_sles_setup_native_conv_diff (int f_id, const char *name, const cs_lnum_t diag_block_size, const cs_lnum_t extra_diag_block_size, const cs_real_t *da, const cs_real_t *xa, bool conv_diff) |
Call sparse linear equation solver setup for convection-diffusion systems. More... | |
cs_sles_convergence_state_t | cs_sles_solve_ccc_fv (cs_sles_t *sc, cs_matrix_t *a, double precision, double r_norm, int *n_iter, double *residual, const cs_real_t *rhs, cs_real_t *vx) |
Call sparse linear equation solver for general colocated cell-centered finite volume scheme. More... | |
cs_sles_convergence_state_t | cs_sles_solve_native (int f_id, const char *name, bool symmetric, cs_lnum_t diag_block_size, cs_lnum_t extra_diag_block_size, const cs_real_t *da, const cs_real_t *xa, double precision, double r_norm, int *n_iter, double *residual, const cs_real_t *rhs, cs_real_t *vx) |
Call sparse linear equation solver using native matrix arrays. More... | |
void | cs_sles_free_native (int f_id, const char *name) |
Free sparse linear equation solver setup using native matrix arrays. More... | |
bool | cs_sles_default_error (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 attempting fallback to alternative solution procedure for sparse linear equation solver. More... | |
void cs_sles_default | ( | int | f_id, |
const char * | name, | ||
const cs_matrix_t * | a | ||
) |
Default definition of a sparse linear equation solver.
[in] | f_id | associated field id, or < 0 |
[in] | name | associated name if f_id < 0, or nullptr |
[in] | a | matrix |
bool cs_sles_default_error | ( | 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 attempting fallback to alternative solution procedure for sparse linear equation solver.
In case of divergence with an iterative solver, this error handler switches to a default preconditioner, then resets the solution vector.
The default error for the solver type handler is then set, in case the solution fails again.
Note that this error handler may rebuild solver contexts, so should not be used in conjunction with shared contexts (such as multigrid ascent/descent contexts), but only for "outer" solvers.
[in,out] | sles | pointer to solver object |
[in] | state | convergence status |
[in] | a | matrix |
[in] | rhs | right hand side |
[in,out] | vx | system solution |
void cs_sles_default_finalize | ( | void | ) |
Default finalization for sparse linear equation solver API.
This includes performance data logging output.
cs_matrix_t * cs_sles_default_get_matrix | ( | int | f_id, |
const char * | name, | ||
cs_lnum_t | db_size, | ||
cs_lnum_t | eb_size, | ||
bool | symmetric | ||
) |
Return pointer to matrix structure matching equation solve.
Some matrix properties (such as assembly options and geometric association) are set immediately, but coefficients are not assigned at this stage.
[in] | f_id | associated field id, or < 0 |
[in] | name | associated name if f_id < 0, or nullptr |
[in] | db_size | block sizes for diagonal |
[in] | eb_size | block sizes for extra diagonal |
[in] | symmetric | indicate if matrix is symmetric |
int cs_sles_default_get_verbosity | ( | int | f_id, |
const char * | name | ||
) |
Return default verbosity associated to a field id, name couple.
[in] | f_id | associated field id, or < 0 |
[in] | name | associated name if f_id < 0, or nullptr |
void cs_sles_default_release_matrix | ( | cs_matrix_t ** | m | ) |
void cs_sles_default_setup | ( | void | ) |
Default setup for sparse linear equation solver API.
This includes setup logging.
void cs_sles_free_native | ( | int | f_id, |
const char * | name | ||
) |
Free sparse linear equation solver setup using native matrix arrays.
[in] | f_id | associated field id, or < 0 |
[in] | name | associated name if f_id < 0, or nullptr |
void cs_sles_setup_native_conv_diff | ( | int | f_id, |
const char * | name, | ||
const cs_lnum_t | diag_block_size, | ||
const cs_lnum_t | extra_diag_block_size, | ||
const cs_real_t * | da, | ||
const cs_real_t * | xa, | ||
bool | conv_diff | ||
) |
Call sparse linear equation solver setup for convection-diffusion systems.
[in] | f_id | associated field id, or < 0 |
[in] | name | associated name if f_id < 0, or nullptr |
[in] | diag_block_size | block sizes for diagonal, or nullptr |
[in] | extra_diag_block_size | block sizes for extra diagonal, or nullptr |
[in] | da | diagonal values (nullptr if zero) |
[in] | xa | extradiagonal values (nullptr if zero) |
[in] | conv_diff | convection-diffusion mode |
cs_sles_convergence_state_t cs_sles_solve_ccc_fv | ( | cs_sles_t * | sc, |
cs_matrix_t * | a, | ||
double | precision, | ||
double | r_norm, | ||
int * | n_iter, | ||
double * | residual, | ||
const cs_real_t * | rhs, | ||
cs_real_t * | vx | ||
) |
Call sparse linear equation solver for general colocated cell-centered finite volume scheme.
The initial solution is assumed to be 0 (and does not need to be initialized before calling this function).
[in] | sc | solver context |
[in] | a | matrix |
[in] | precision | solver precision |
[in] | r_norm | residual normalization |
[out] | n_iter | number of "equivalent" iterations |
[out] | residual | residual |
[in] | rhs | right hand side |
[out] | vx | system solution |
cs_sles_convergence_state_t cs_sles_solve_native | ( | int | f_id, |
const char * | name, | ||
bool | symmetric, | ||
cs_lnum_t | diag_block_size, | ||
cs_lnum_t | extra_diag_block_size, | ||
const cs_real_t * | da, | ||
const cs_real_t * | xa, | ||
double | precision, | ||
double | r_norm, | ||
int * | n_iter, | ||
double * | residual, | ||
const cs_real_t * | rhs, | ||
cs_real_t * | vx | ||
) |
Call sparse linear equation solver using native matrix arrays.
[in] | f_id | associated field id, or < 0 |
[in] | name | associated name if f_id < 0, or nullptr |
[in] | symmetric | indicates if matrix coefficients are symmetric |
[in] | diag_block_size | block sizes for diagonal |
[in] | extra_diag_block_size | block sizes for extra diagonal |
[in] | da | diagonal values (nullptr if zero) |
[in] | xa | extradiagonal values (nullptr if zero) |
[in] | precision | solver precision |
[in] | r_norm | residual normalization |
[out] | n_iter | number of "equivalent" iterations |
[out] | residual | residual |
[in] | rhs | right hand side |
[in,out] | vx | system solution |