strumpack::SPOptions< scalar_t > Class Template Reference

Options for the sparse solver. More...

#include <StrumpackOptions.hpp>

Public Types

using real_t = typename RealType< scalar_t >::value_type
 

Public Member Functions

 SPOptions ()
 
 SPOptions (int argc, const char *const argv[])
 
void set_verbose (bool verbose)
 
void set_maxit (int maxit)
 
void set_rel_tol (real_t rtol)
 
void set_abs_tol (real_t atol)
 
void set_Krylov_solver (KrylovSolver s)
 
void set_gmres_restart (int m)
 
void set_GramSchmidt_type (GramSchmidtType t)
 
void set_reordering_method (ReorderingStrategy m)
 
void set_nd_param (int nd_param)
 
void set_nd_planar_levels (int nd_planar_levels)
 
void set_dimensions (int nx, int ny=1, int nz=1)
 
void set_nx (int nx)
 
void set_ny (int ny)
 
void set_nz (int nz)
 
void set_components (int components)
 
void set_separator_width (int width)
 
void enable_METIS_NodeNDP ()
 
void disable_METIS_NodeNDP ()
 
void enable_METIS_NodeND ()
 
void disable_METIS_NodeND ()
 
void enable_MUMPS_SYMQAMD ()
 
void disable_MUMPS_SYMQAMD ()
 
void enable_agg_amalg ()
 
void disable_agg_amalg ()
 
void set_matching (MatchingJob job)
 
void enable_assembly_tree_log ()
 
void disable_assembly_tree_log ()
 
void set_compression (CompressionType c)
 
void set_compression_rel_tol (real_t rtol)
 
void set_compression_abs_tol (real_t atol)
 
void set_compression_min_sep_size (int s)
 
void set_hss_min_sep_size (int s)
 
void set_hodlr_min_sep_size (int s)
 
void set_blr_min_sep_size (int s)
 
void set_lossy_min_sep_size (int s)
 
void set_compression_min_front_size (int s)
 
void set_hss_min_front_size (int s)
 
void set_hodlr_min_front_size (int s)
 
void set_blr_min_front_size (int s)
 
void set_lossy_min_front_size (int s)
 
void set_compression_leaf_size (int s)
 
void set_separator_ordering_level (int l)
 
void enable_replace_tiny_pivots ()
 
void disable_replace_tiny_pivots ()
 
void set_pivot_threshold (real_t thresh)
 
void set_write_root_front (bool b)
 
void enable_gpu ()
 
void disable_gpu ()
 
void set_gpu_streams (int s)
 
void enable_openmp_tree ()
 
void disable_openmp_tree ()
 
void set_lossy_precision (int p)
 
void set_print_compressed_front_stats (bool b)
 
void set_proportional_mapping (ProportionalMapping pmap)
 
bool verbose () const
 
int maxit () const
 
real_t rel_tol () const
 
real_t abs_tol () const
 
KrylovSolver Krylov_solver () const
 
int gmres_restart () const
 
GramSchmidtType GramSchmidt_type () const
 
ReorderingStrategy reordering_method () const
 
int nd_param () const
 
int nd_planar_levels () const
 
int nx () const
 
int ny () const
 
int nz () const
 
int components () const
 
int separator_width () const
 
bool use_METIS_NodeNDP () const
 
bool use_METIS_NodeND () const
 
bool use_MUMPS_SYMQAMD () const
 
bool use_agg_amalg () const
 
MatchingJob matching () const
 
bool log_assembly_tree () const
 
CompressionType compression () const
 
real_t compression_rel_tol (int l=0) const
 
real_t compression_abs_tol (int l=0) const
 
int compression_min_sep_size (int l=0) const
 
int hss_min_sep_size () const
 
int hodlr_min_sep_size () const
 
int blr_min_sep_size () const
 
int lossy_min_sep_size () const
 
int compression_min_front_size (int l=0) const
 
int hss_min_front_size () const
 
int hodlr_min_front_size () const
 
int blr_min_front_size () const
 
int lossy_min_front_size () const
 
int compression_leaf_size (int l=0) const
 
int separator_ordering_level () const
 
bool indirect_sampling () const
 
bool replace_tiny_pivots () const
 
real_t pivot_threshold () const
 
bool write_root_front () const
 
bool use_gpu () const
 
bool use_openmp_tree () const
 
int gpu_streams () const
 
int lossy_precision () const
 
bool print_compressed_front_stats () const
 
ProportionalMapping proportional_mapping () const
 
const HSS::HSSOptions< scalar_t > & HSS_options () const
 
HSS::HSSOptions< scalar_t > & HSS_options ()
 
const BLR::BLROptions< scalar_t > & BLR_options () const
 
BLR::BLROptions< scalar_t > & BLR_options ()
 
const HODLR::HODLROptions< scalar_t > & HODLR_options () const
 
HODLR::HODLROptions< scalar_t > & HODLR_options ()
 
void set_from_command_line ()
 
void set_from_command_line (int argc, const char *const *cargv)
 
void describe_options () const
 

Detailed Description

template<typename scalar_t>
class strumpack::SPOptions< scalar_t >

Options for the sparse solver.

This sparse solver object also stores an object with HSS options (HSS_options), one with BLR options (BLR_options) and one with HODLR options (HODLR_options), since HSS, BLR and HODLR compression can be used in the sparse solver.

Running with -h or –help will print a list of options when the set_from_command_line() routine is called.

Template Parameters
scalar_tcan be float, double, std::complex<float> or std::complex<double>, should be the same as used for the StrumpackSparseSolver object

Member Typedef Documentation

◆ real_t

template<typename scalar_t >
using strumpack::SPOptions< scalar_t >::real_t = typename RealType<scalar_t>::value_type

real_t is the real type corresponding to the (possibly complex) scalar_t template parameter

Constructor & Destructor Documentation

◆ SPOptions() [1/2]

template<typename scalar_t >
strumpack::SPOptions< scalar_t >::SPOptions ( )
inline

Default constructor, initializing all options to their default values.

◆ SPOptions() [2/2]

template<typename scalar_t >
strumpack::SPOptions< scalar_t >::SPOptions ( int  argc,
const char *const  argv[] 
)
inline

Constructor, initializing all options to their default values. This will store a copy of *argv[]. It will not yet parse the options!! The options are only parsed when calling set_from_command_line(). This allows the user of this class to set certain options using the member functions of this class, and then overwrite those with command line arguments later.

Parameters
argcnumber of arguments (number of char* in argv)
argvinout from the command line with list of options

Member Function Documentation

◆ abs_tol()

template<typename scalar_t >
real_t strumpack::SPOptions< scalar_t >::abs_tol ( ) const
inline

Get the absolute tolerance to be used in the iterative solver.

See also
set_abs_tol()

◆ BLR_options() [1/2]

template<typename scalar_t >
BLR::BLROptions< scalar_t > & strumpack::SPOptions< scalar_t >::BLR_options ( )
inline

Get a reference to an object holding various options pertaining to the BLR code, and data structures.

◆ BLR_options() [2/2]

template<typename scalar_t >
const BLR::BLROptions< scalar_t > & strumpack::SPOptions< scalar_t >::BLR_options ( ) const
inline

Get a (const) reference to an object holding various options pertaining to the BLR code, and data structures.

◆ components()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::components ( ) const
inline

Get the currently specified number of components (DoF's) per mesh point.

See also
set_components()

◆ compression()

template<typename scalar_t >
CompressionType strumpack::SPOptions< scalar_t >::compression ( ) const
inline

Get the type of compression to use.

◆ compression_abs_tol()

template<typename scalar_t >
real_t strumpack::SPOptions< scalar_t >::compression_abs_tol ( int  l = 0) const
inline

Return the absolute compression tolerance used for the currently selected low-rank method, either HODLR, HSS or BLR. If NONE, LOSSY or LOSSLESS compression are selected, then this returns 0.

See also
set_compression_rel_tol, set_lossy_compression

◆ compression_leaf_size()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::compression_leaf_size ( int  l = 0) const
inline

Get the leaf size used in the rank-structured format used for compression. This will depend on which type of compression is selected.

See also
set_compression(), set_compression_leaf_size(), compression_min_sep_size()

◆ compression_min_front_size()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::compression_min_front_size ( int  l = 0) const
inline

Get the minimum size of a front to enable compression. This will depend on which type of compression is selected.

See also
set_compression(), set_compression_min_sep_size(), compression_min_front_size()

◆ compression_min_sep_size()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::compression_min_sep_size ( int  l = 0) const
inline

Get the minimum size of a separator to enable compression. This will depend on which type of compression is selected.

See also
set_compression(), set_compression_min_sep_size(), compression_min_front_size()

◆ compression_rel_tol()

template<typename scalar_t >
real_t strumpack::SPOptions< scalar_t >::compression_rel_tol ( int  l = 0) const
inline

Return the relative compression tolerance used for the currently selected low-rank method, either HODLR, HSS or BLR. If NONE, LOSSY or LOSSLESS compression are selected, then this returns 0.

See also
set_compression_rel_tol, set_lossy_compression

◆ describe_options()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::describe_options ( ) const

Print an overview of all supported options. Not including any HSS/BLR specific options.

◆ disable_agg_amalg()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::disable_agg_amalg ( )
inline

Disbale aggressive amalgamation of nodes into supernodes inside MUMPS_SYMQMAD. This is only relevant when MUMPS_SYMQAMD is enabled.

See also
enable_agg_amalg(), enable_MUMPS_SYMQAMD(), disable_MUMPS_SYMQAMD()

◆ disable_assembly_tree_log()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::disable_assembly_tree_log ( )
inline

Do not log the assembly tree to a file. Logging of the tree is currently not supported.

◆ disable_gpu()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::disable_gpu ( )
inline

Disable GPU off-loading.

◆ disable_METIS_NodeND()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::disable_METIS_NodeND ( )
inline

Do not use the routine METIS_NodeND, but instead use the undocumented routine METIS_NodeNDP.

See also
enable_METIS_NodeNDP(), disable_METIS_NodeNDP(), enable_METIS_NodeND()

◆ disable_METIS_NodeNDP()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::disable_METIS_NodeNDP ( )
inline

Disable use of the routine METIS_NodeNDP, and instead use METIS_NodeND. METIS_NodeNDP is a non-documented Metis routine to perform nested dissection which (unlike METIS_NodeND) also return the separator tree.

See also
enable_METIS_NodeNDP(), enable_METIS_NodeND(), disable_METIS_NodeND()

◆ disable_MUMPS_SYMQAMD()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::disable_MUMPS_SYMQAMD ( )
inline

Do not use the SYMQAMD routine (provided by the MUMPS folks) to construct the supernodal tree from the elimination tree. In some cases SYMQAMD can reduce the amount of fill. Use STRUMPACK's own routine instead.

See also
enable_MUMPS_SYMQAMD(), enable_agg_amalg(), disable_agg_amalg()

◆ disable_openmp_tree()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::disable_openmp_tree ( )
inline

Disable OpenMP tasking traversal of the supernodal tree in the sparse solver. This reduces the peak memory requirements, but scales poorer with OpenMP threads.

◆ disable_replace_tiny_pivots()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::disable_replace_tiny_pivots ( )
inline

Disable replacement of tiny pivots.

See also
enable_replace_tiny_pivots()

◆ enable_agg_amalg()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::enable_agg_amalg ( )
inline

When using MUMPS_SYMQAMD, enable aggressive amalgamation of nodes into supernodes. This is only used when MUMPS_SYMQAMD is enabled.

See also
disable_agg_amalg(), enable_MUMPS_SYMQAMD(), disable_MUMPS_SYMQAMD()

◆ enable_assembly_tree_log()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::enable_assembly_tree_log ( )
inline

Log the assembly tree to a file. Currently not supported.

◆ enable_gpu()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::enable_gpu ( )
inline

Enable off-loading to the GPU. This only works when STRUMPACK was configured with GPU support (through CUDA, MAGMA or SLATE)

◆ enable_METIS_NodeND()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::enable_METIS_NodeND ( )
inline

Use the routine METIS_NodeND instead of the undocumented routine METIS_NodeNDP.

See also
enable_METIS_NodeNDP(), disable_METIS_NodeNDP(), disable_METIS_NodeND()

◆ enable_METIS_NodeNDP()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::enable_METIS_NodeNDP ( )
inline

Enable use of the routine METIS_NodeNDP, instead of METIS_NodeND. METIS_NodeNDP is a non-documented Metis routine to perform nested dissection which (unlike METIS_NodeND) also return the separator tree.

See also
disable_METIS_NodeNDP(), enable_METIS_NodeND(), disable_METIS_NodeND()

◆ enable_MUMPS_SYMQAMD()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::enable_MUMPS_SYMQAMD ( )
inline

Use the SYMQAMD routine (provided by the MUMPS folks) to construct the supernodal tree from the elimination tree. In some cases this can reduce the amount of fill.

See also
disable_MUMPS_SYMQAMD(), enable_agg_amalg(), disable_agg_amalg()

◆ enable_openmp_tree()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::enable_openmp_tree ( )
inline

Enable OpenMP tasking traversal of the supernodal tree in the sparse solver. This requires more (peak) memory, but scales better with OpenMP threads.

◆ enable_replace_tiny_pivots()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::enable_replace_tiny_pivots ( )
inline

Enable replacing of small pivot values with a larger value. This can prevent to numerical factorization to fail completely, but the resulting solve might be inaccurate and hence require more iterations. If the outer iterative solver fails to converge, you can switch from iterative refinement (which is used with the AUTO solver strategy if compression is disabled), to PGMRES (preconditioned!).

If you encounter numerical issues during the factorization (such as small pivots, failure in LU factorization), you can also try a different matching (MatchingJob::MAX_DIAGONAL_PRODUCT_SCALING is the most robust option).

See also
disable_replace_tiny_pivots(), set_matching()

◆ gmres_restart()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::gmres_restart ( ) const
inline

Get the GMRES restart length.

See also
set_gmres_restart()

◆ gpu_streams()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::gpu_streams ( ) const
inline

Returns the number of GPU streams to use.

◆ GramSchmidt_type()

template<typename scalar_t >
GramSchmidtType strumpack::SPOptions< scalar_t >::GramSchmidt_type ( ) const
inline

Get the Gram-Schmidth orthogonalization type used in GMRES.

See also
set_GramSchmidth_type()

◆ HODLR_options() [1/2]

template<typename scalar_t >
HODLR::HODLROptions< scalar_t > & strumpack::SPOptions< scalar_t >::HODLR_options ( )
inline

Get a reference to an object holding various options pertaining to the HODLR code, and data structures.

◆ HODLR_options() [2/2]

template<typename scalar_t >
const HODLR::HODLROptions< scalar_t > & strumpack::SPOptions< scalar_t >::HODLR_options ( ) const
inline

Get a (const) reference to an object holding various options pertaining to the HODLR code, and data structures.

◆ HSS_options() [1/2]

template<typename scalar_t >
HSS::HSSOptions< scalar_t > & strumpack::SPOptions< scalar_t >::HSS_options ( )
inline

Get a reference to an object holding various options pertaining to the HSS code, and data structures.

◆ HSS_options() [2/2]

template<typename scalar_t >
const HSS::HSSOptions< scalar_t > & strumpack::SPOptions< scalar_t >::HSS_options ( ) const
inline

Get a (const) reference to an object holding various options pertaining to the HSS code, and data structures.

◆ indirect_sampling()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::indirect_sampling ( ) const
inline

Is indirect sampling for HSS construction enabled?

◆ Krylov_solver()

template<typename scalar_t >
KrylovSolver strumpack::SPOptions< scalar_t >::Krylov_solver ( ) const
inline

Get the type of iterative solver to be used as outer solver.

See also
set_Krylov_solver()

◆ log_assembly_tree()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::log_assembly_tree ( ) const
inline

Should we log the assembly tree? Currently not supported.

◆ lossy_precision()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::lossy_precision ( ) const
inline

Returns the precision for lossy compression.

◆ matching()

template<typename scalar_t >
MatchingJob strumpack::SPOptions< scalar_t >::matching ( ) const
inline

Get the matching job to use for numerical stability reordering.

See also
set_matching()

◆ maxit()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::maxit ( ) const
inline

Get the maximum number of allowed iterative solver iterations.

See also
set_maxit()

◆ nd_param()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::nd_param ( ) const
inline

Return the nested-dissection recursion ending parameter.

See also
set_nd_param()

◆ nd_planar_levels()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::nd_planar_levels ( ) const
inline

Return the number of levels in nested-dissection for which the separators are parallel instead of split perpedicular to the longest direction. This is to reduce the ranks in the F12/F21 blocks when using compression.

◆ nx()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::nx ( ) const
inline

Get the specified nx mesh dimension.

See also
set_nx()

◆ ny()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::ny ( ) const
inline

Get the specified ny mesh dimension.

See also
set_ny()

◆ nz()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::nz ( ) const
inline

Get the specified nz mesh dimension.

See also
set_nz()

◆ pivot_threshold()

template<typename scalar_t >
real_t strumpack::SPOptions< scalar_t >::pivot_threshold ( ) const
inline

Get the minimum pivot value. If the option replace_tiny_pivots is enabled (using enable_replace_tiny_pivots()), the all pivots smaller than this threshold value will be replaced by this threshold value.

See also
enable_replace_tiny_pivots() set_pivot_threshold()

◆ print_compressed_front_stats()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::print_compressed_front_stats ( ) const
inline

Info about the stats of the root front will be printed to std::cout

◆ proportional_mapping()

template<typename scalar_t >
ProportionalMapping strumpack::SPOptions< scalar_t >::proportional_mapping ( ) const
inline

Get the type of proportional mapping to be used.

◆ rel_tol()

template<typename scalar_t >
real_t strumpack::SPOptions< scalar_t >::rel_tol ( ) const
inline

Get the relative tolerance to be used in the iterative solver.

See also
set_rel_tol()

◆ reordering_method()

template<typename scalar_t >
ReorderingStrategy strumpack::SPOptions< scalar_t >::reordering_method ( ) const
inline

Get the currently set fill reducing reordering method.

See also
set_reordering_method()

◆ replace_tiny_pivots()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::replace_tiny_pivots ( ) const
inline

Check whether replacement of tiny pivots is enabled.

See also
enable_replace_tiny_pivots()

◆ separator_ordering_level()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::separator_ordering_level ( ) const
inline

Get the current value of the number of additional links to include in the separator before reordering.

See also
set_separator_ordering_level()

◆ separator_width()

template<typename scalar_t >
int strumpack::SPOptions< scalar_t >::separator_width ( ) const
inline

Get the currently specified width of a separator.

See also
set_separator_width()

◆ set_abs_tol()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_abs_tol ( real_t  atol)
inline

Set the absolute tolerance to use in the iterative solvers.

Parameters
abs_tolabsolute tolerance

◆ set_components()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_components ( int  components)
inline

Set the number of components per gridpoint. This is only useful when the sparse matrix was generated by a stencil on a regular 1d, 2d or 3d mesh. The degrees of freedom for a single gridpoint should be ordered consecutively.

Parameters
componentsnumber of components per gridpoint
See also
set_dimensions(), set_separator_width(), ReorderingStrategy::GEOMETRIC, set_reordering_method()

◆ set_compression()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_compression ( CompressionType  c)
inline

Set the type of rank-structured compression to use.

Parameters
ccompression type
See also
set_compression_min_sep_size(), set_compression_min_front_size()

◆ set_compression_abs_tol()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_compression_abs_tol ( real_t  atol)
inline

Set the absolute compression tolerance to be used for low-rank compression. This currently affects BLR, HSS, HODLR, HODBF, Butterfly. It does not affect the lossy compression, see set_lossy_precision. The same tolerances can also be set individually by calling set_rel_tol on the returned objects from either HODLR_options(), BLR_options() or HSS_options().

Parameters
rtolthe relative low-rank compression tolerance
See also
set_compression_rel_tol, set_lossy_precision, compression_abs_tol

◆ set_compression_leaf_size()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_compression_leaf_size ( int  s)
inline

Set the leaf size used by any of the rank-structured formats.

See also
HSS_options(), BLR_options(), HODLR_options()

◆ set_compression_min_front_size()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_compression_min_front_size ( int  s)
inline

Set the minimum size of frontal matrices (dense submatrices of the sparse triangular factors), for which to use compression (when HSS/BLR/HODLR/LOSSY compression has been enabled by the user. Most fronts will be quite small, and will not have benefit from compression.

See also
set_compression(), set_compression_min_front_size()

◆ set_compression_min_sep_size()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_compression_min_sep_size ( int  s)
inline

Set the minimum size of the top left part of frontal matrices (dense submatrices of the sparse triangular factors), ie, the part corresponding to the separators, for which to use compression (when HSS/BLR/HODLR/LOSSY compression has been enabled by the user. Most fronts will be quite small, and will not have benefit from compression.

See also
set_compression(), set_compression_min_front_size()

◆ set_compression_rel_tol()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_compression_rel_tol ( real_t  rtol)
inline

Set the relative compression tolerance to be used for low-rank compression. This currently affects BLR, HSS, HODLR, HODBF, Butterfly. It does not affect the lossy compression, see set_lossy_precision. The same tolerances can also be set individually by calling set_rel_tol on the returned objects from either HODLR_options(), BLR_options() or HSS_options().

Parameters
rtolthe relative low-rank compression tolerance
See also
set_compression_abs_tol, set_lossy_precision, compression_rel_tol

◆ set_dimensions()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_dimensions ( int  nx,
int  ny = 1,
int  nz = 1 
)
inline

Set the mesh dimensions. This is only useful when the sparse matrix was generated by a stencil on a regular 1d, 2d or 3d mesh. The stencil can be multiple gridpoints wide, and there can be multiple degrees of freedom per gridpoint. The degrees of freedom for a single gridpoint should be ordered consecutively. When the dimensions of the grid are specified, along with the width of the stencil and the number of degrees of freedom per gridpoint, one can use the GEOMETRIC option for the fill-reducing reordering.

Parameters
nxdimension of the grid along the first axis
nydimension of the grid along the second axis, should not be specified for a 1d mesh
nzdimension of the grid along the third axis, should not be specified for a 2d mesh
See also
set_components(), set_separator_width(), ReorderingStrategy::GEOMETRIC, set_reordering_method()

◆ set_from_command_line() [1/2]

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_from_command_line ( )
inline

Parse the command line options that were passed to this object in the constructor. Run the code with -h or –help and call this routine to see a list of supported options.

◆ set_from_command_line() [2/2]

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_from_command_line ( int  argc,
const char *const *  cargv 
)

Parse command line options. These options will not be modified. The options can also contain HSS, BLR or HODLR specific options, they will be parsed by the HSS::HSSOptions and BLR::BLROptions, HODLR::HODLROptions objects returned by HSS_options(), BLR_options() and HODLR_options() respectively. Run the code with -h or –help and call this routine to see a list of supported options.

Parameters
argcnumber of arguments in the argv array
argvlist of options

◆ set_gmres_restart()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_gmres_restart ( int  m)
inline

Set the GMRES restart length

Parameters
mGMRES restart, should be > 0 and not too crazy big

◆ set_gpu_streams()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_gpu_streams ( int  s)
inline

Set the number of (CUDA) streams to be used in the code.

◆ set_GramSchmidt_type()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_GramSchmidt_type ( GramSchmidtType  t)
inline

Set the type of Gram-Schmidt orthogonalization to use in GMRES

Parameters
tGram-Schmidt type to use in GMRES

◆ set_Krylov_solver()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_Krylov_solver ( KrylovSolver  s)
inline

Select a Krylov outer solver. Note that the versions GMRES and BICGSTAB are not preconditioned! Use PREC_GMRES and PREC_BICGSTAB instead, as they will use STRUMPACK's sparse solver (possibly with compression) as a preconditioner. Setting this to DIRECT will not use any iterative solver, and instead only perform a solve with the (incomplete/compressed) LU factorization.

Parameters
souter, iterative solver to use
See also
set_compression()

◆ set_lossy_precision()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_lossy_precision ( int  p)
inline

Set the precision for lossy compression.

◆ set_matching()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_matching ( MatchingJob  job)
inline

Specify the job type for the column ordering for stability. This ordering is computed using a maximum matching algorithm, to try to get nonzero (as large as possible) values on the diagonal to improve numerica stability. Possibly, this can also scale the rows and columns of the matrix.

Parameters
jobtype of matching to perform

◆ set_maxit()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_maxit ( int  maxit)
inline

Set the maximum number of iterations to use in any of the iterative solvers.

◆ set_nd_param()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_nd_param ( int  nd_param)
inline

Set the parameter used in nested dissection to determine when to stop the nested-dissection recursion. Larger values may lead to less nodes in the separator tree, which eliminates some overhead, but typically leads to more fill.

◆ set_nd_planar_levels()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_nd_planar_levels ( int  nd_planar_levels)
inline

Set the number of levels in nested-dissection for which the separators are parallel instead of perpedicular to the longest direction. This is to reduce the ranks in the F12/F21 blocks when using compression.

◆ set_nx()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_nx ( int  nx)
inline

Set the mesh dimension along the x-axis. This is only useful when the sparse matrix was generated by a stencil on a regular 1d, 2d or 3d mesh.

See also
set_dimensions(), set_ny(), set_nz(), set_components(), set_separator_width(), ReorderingStrategy::GEOMETRIC, set_reordering_method()

◆ set_ny()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_ny ( int  ny)
inline

Set the mesh dimension along the y-axis. This is only useful when the sparse matrix was generated by a stencil on a regular 1d, 2d or 3d mesh. For a 1d mesh this should not be specified, defaults to 1.

See also
set_dimensions(), set_nx(), set_nz(), set_components(), set_separator_width(), ReorderingStrategy::GEOMETRIC, set_reordering_method()

◆ set_nz()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_nz ( int  nz)
inline

Set the mesh dimension along the z-axis. This is only useful when the sparse matrix was generated by a stencil on a regular 1d, 2d or 3d mesh. For a 2d mesh this should not be specified, defaults to 1.

See also
set_dimensions(), set_nx(), set_ny(), set_components(), set_separator_width(), ReorderingStrategy::GEOMETRIC, set_reordering_method()

◆ set_pivot_threshold()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_pivot_threshold ( real_t  thresh)
inline

Set the minimum pivot value. If the option replace_tiny_pivots is enabled (using enable_replace_tiny_pivots()), the all pivots smaller than this threshold value will be replaced by this threshold value. This should not be set by the user directly, but is set in the solver.

See also
enable_replace_tiny_pivots()

◆ set_print_compressed_front_stats()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_print_compressed_front_stats ( bool  b)
inline

Print statistics, about ranks, memory etc, for the root front only.

◆ set_proportional_mapping()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_proportional_mapping ( ProportionalMapping  pmap)
inline

Set the type of proportional mapping.

◆ set_rel_tol()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_rel_tol ( real_t  rtol)
inline

Set the relative tolerance to use in the iterative solvers. This will typically be the relative residual decrease.

Parameters
rel_tolrelative tolerance

◆ set_reordering_method()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_reordering_method ( ReorderingStrategy  m)
inline

Set the sparse fill-reducing reordering. This can greatly affect the memory usage and factorization time. However, note that most reordering routines are provided by third party libraries, such as Metis and Scotch. In order to use those, STRUMPACK needs to be configured with support for those libraries. Some reorderings only work when using the distributed memory solvers (StrumpackSparseSolverMPI or StrumpackSparseSolverMPIDist).

Parameters
mfill reducing reordering

◆ set_separator_ordering_level()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_separator_ordering_level ( int  l)
inline

This is used in the reordering of the separators, which is useful when using compression (HSS, BLR) to reduce the ranks and to define the block-sizes in the rank structured matrix format. A value of 1 means that additional edges are introduced between nodes in the separators, edges corresponding to 2-step connections in the graph of the original matrix, before passing the graph to a graph partitioner. These additional edges are useful beause sometimes the separator graph can become unconnected (although the points are connected in the original graph). A value of 0 means not to add these additional edges.

You should probably leave this at 1, unless separator reordering takes a long time, in which case you can try 0.

◆ set_separator_width()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_separator_width ( int  width)
inline

Set the width of the separator, which depends on the width of the stencil. For instance for a 1d 3-point stencil, the separator is just a single point, for a 2d 5-point stencil the separator is a single line and for a 3d stencil, the separator is a single plane, and hence, for all these cases, the separator width is 1. Wider stencils need a wider separator, for instance a stencil with 5 points in each dimension needs a separator width of 3. This is only useful when the sparse matrix was generated by a stencil on a regular 1d, 2d or 3d mesh. The degrees of freedom for a single gridpoint should be ordered consecutively.

Parameters
widthwidth of the separator
See also
set_dimensions(), set_components(), ReorderingStrategy::GEOMETRIC, set_reordering_method()

◆ set_verbose()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_verbose ( bool  verbose)
inline

Set verbose to true/false, ie, allow the sparse solver to print out progress information, statistics on time, flops, memory usage etc. to cout. Warnings will go to cerr regardless of the verbose options set here.

◆ set_write_root_front()

template<typename scalar_t >
void strumpack::SPOptions< scalar_t >::set_write_root_front ( bool  b)
inline

Dump the root front to a set of files, one for each rank. This will only have affect when running with more than one MPI rank, and without compression.

◆ use_agg_amalg()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::use_agg_amalg ( ) const
inline

Is aggressive amalgamation enabled? (only used when MUMPS_SYMQAMD is enabled)

See also
enable_agg_amalg(), enable_MUMPS_SYMQAMD()

◆ use_gpu()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::use_gpu ( ) const
inline

Check wheter or not to use GPU off-loading.

◆ use_METIS_NodeND()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::use_METIS_NodeND ( ) const
inline

Is use of METIS_NodeND enabled? (instead of METIS_NodeNDP)

See also
enable_METIS_NodeND()

◆ use_METIS_NodeNDP()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::use_METIS_NodeNDP ( ) const
inline

Is use of METIS_NodeNDP enabled? (instead of METIS_NodeND)

See also
enable_METIS_NodeNDP()

◆ use_MUMPS_SYMQAMD()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::use_MUMPS_SYMQAMD ( ) const
inline

Is MUMPS_SYMQAMD enabled?

See also
enable_MUMPS_SYMQAMD()

◆ use_openmp_tree()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::use_openmp_tree ( ) const
inline

Check wheter or not to use OpenMP tree traversal is the sparse solver.

◆ verbose()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::verbose ( ) const
inline

Check if verbose output is enabled.

See also
set_verbose()

◆ write_root_front()

template<typename scalar_t >
bool strumpack::SPOptions< scalar_t >::write_root_front ( ) const
inline

The root front will be written to a file.


The documentation for this class was generated from the following file: