My Project
programmer's documentation
|
#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_blas.h"
#include "cs_file.h"
#include "cs_log.h"
#include "cs_halo.h"
#include "cs_mesh.h"
#include "cs_matrix.h"
#include "cs_matrix_default.h"
#include "cs_matrix_util.h"
#include "cs_post.h"
#include "cs_timer.h"
#include "cs_time_plot.h"
#include "cs_sles.h"
#include "cs_sles_pc.h"
#include "cs_sles_it.h"
Functions | |
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. More... | |
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. 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_log (const void *context, cs_log_t log_type) |
Log sparse linear equation solver info. 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, cs_halo_rotation_t rotation_mode, 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... | |
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_t * | cs_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... | |
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... | |
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, cs_halo_rotation_t rotation_mode, 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... | |
void | cs_sles_it_assign_plot (cs_sles_it_t *context, cs_time_plot_t *time_plot, int time_stamp) |
Assign existing time plot to iterative sparse linear equation solver. More... | |
Iterative linear solvers
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.
[in,out] | context | pointer to iterative solver info and context |
[in,out] | order | pointer to ordering array |
void cs_sles_it_assign_plot | ( | cs_sles_it_t * | context, |
cs_time_plot_t * | time_plot, | ||
int | time_stamp | ||
) |
Assign existing time plot to iterative sparse linear equation solver.
This is useful mainly when a time plot has a longer lifecycle than the linear solver context, such as inside a multigrid solver.
[in,out] | context | pointer to iterative solver info and context |
[in] | time_plot | pointer to time plot structure |
[in] | time_stamp | associated time stamp |
void* cs_sles_it_copy | ( | const void * | context | ) |
Create iterative sparse linear system solver info and context based on existing info and context.
[in] | context | pointer to reference info and context (actual type: cs_sles_it_t *) |
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.
[in] | solver_type | type of solver (PCG, Jacobi, ...) |
[in] | poly_degree | preconditioning polynomial degree (0: diagonal; -1: non-preconditioned; see Iterative linear solvers. for details) |
[in] | n_max_iter | maximum number of iterations |
[in] | update_stats | automatic solver statistics indicator |
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.
[in] | f_id | associated field id, or < 0 |
[in] | name | associated name if f_id < 0, or NULL |
[in] | solver_type | type of solver (PCG, Jacobi, ...) |
[in] | poly_degree | preconditioning polynomial degree (0: diagonal; -1: non-preconditioned) |
[in] | n_max_iter | maximum number of iterations |
void cs_sles_it_destroy | ( | void ** | context | ) |
Destroy iterative sparse linear system solver info and context.
[in,out] | context | pointer to iterative solver info and context (actual type: cs_sles_it_t **) |
bool cs_sles_it_error_post_and_abort | ( | cs_sles_t * | sles, |
cs_sles_convergence_state_t | state, | ||
const cs_matrix_t * | a, | ||
cs_halo_rotation_t | rotation_mode, | ||
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.
[in,out] | sles | pointer to solver object |
[in] | state | convergence state |
[in] | a | matrix |
[in] | rotation_mode | halo update option for rotational periodicity |
[in] | rhs | right hand side |
[in,out] | vx | system solution |
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.
[in,out] | context | pointer to iterative solver info and context (actual type: cs_sles_it_t *) |
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.
[in] | context | pointer to iterative solver info and context |
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.
[in] | context | pointer to iterative solver info and context |
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.
cs_sles_it_type_t cs_sles_it_get_type | ( | const cs_sles_it_t * | context | ) |
Return iterative solver type.
[in] | context | pointer to iterative solver info and context |
void cs_sles_it_log | ( | const void * | context, |
cs_log_t | log_type | ||
) |
Log sparse linear equation solver info.
[in] | context | pointer to iterative solver info and context (actual type: cs_sles_it_t *) |
[in] | log_type | log type |
void cs_sles_it_log_parallel_options | ( | void | ) |
Log the current global settings relative to parallelism.
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.
[in,out] | context | pointer to iterative solver info and context |
[in] | threshold | convergence level under which fallback is used |
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.
[in] | threshold | mean number of rows per active rank under which the single-reduction variant will be used |
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.
[in,out] | context | pointer to iterative solver info and context |
[in] | base_name | base plot name to activate, NULL otherwise |
[in] | use_iteration | if true, use iteration as time stamp otherwise, use wall clock time |
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).
[in,out] | context | pointer to iterative solver info and context |
[in] | shareable | pointer to iterative solver info and context whose context may be shared |
void cs_sles_it_setup | ( | void * | context, |
const char * | name, | ||
const cs_matrix_t * | a, | ||
int | verbosity | ||
) |
Setup iterative sparse linear equation solver.
[in,out] | context | pointer to iterative solver info and context (actual type: cs_sles_it_t *) |
[in] | name | pointer to system name |
[in] | a | associated matrix |
[in] | verbosity | associated verbosity |
cs_sles_convergence_state_t cs_sles_it_solve | ( | void * | context, |
const char * | name, | ||
const cs_matrix_t * | a, | ||
int | verbosity, | ||
cs_halo_rotation_t | rotation_mode, | ||
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.
[in,out] | context | pointer to iterative solver info and context (actual type: cs_sles_it_t *) |
[in] | name | pointer to system name |
[in] | a | matrix |
[in] | verbosity | associated verbosity |
[in] | rotation_mode | halo update option for rotational periodicity |
[in] | precision | solver precision |
[in] | r_norm | residue normalization |
[out] | n_iter | number of "equivalent" iterations |
[out] | residue | residue |
[in] | rhs | right hand side |
[in,out] | vx | system solution |
[in] | aux_size | number of elements in aux_vectors (in bytes) |
aux_vectors | optional working area (internal allocation if NULL) |
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.
[in] | src | pointer to source info and context |
[in,out] | dest | pointer to destination info and context |
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.
[in,out] | context | pointer to iterative solver info and context |
[in,out] | pc | pointer to preconditoner context |