Options for the sparse solver. More...
#include <StrumpackOptions.hpp>
Public Types | |
using | real_t = typename RealType< scalar_t >::value_type |
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.
scalar_t | can be float, double, std::complex<float> or std::complex<double>, should be the same as used for the StrumpackSparseSolver object |
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
|
inline |
Default constructor, initializing all options to their default values.
|
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.
argc | number of arguments (number of char* in argv) |
argv | inout from the command line with list of options |
|
inline |
Get the absolute tolerance to be used in the iterative solver.
|
inline |
Get a reference to an object holding various options pertaining to the BLR code, and data structures.
|
inline |
Get a (const) reference to an object holding various options pertaining to the BLR code, and data structures.
|
inline |
Get the currently specified number of components (DoF's) per mesh point.
|
inline |
Get the type of compression to use.
|
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.
|
inline |
Get the leaf size used in the rank-structured format used for compression. This will depend on which type of compression is selected.
|
inline |
Get the minimum size of a front to enable compression. This will depend on which type of compression is selected.
|
inline |
Get the minimum size of a separator to enable compression. This will depend on which type of compression is selected.
|
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.
void strumpack::SPOptions< scalar_t >::describe_options | ( | ) | const |
Print an overview of all supported options. Not including any HSS/BLR specific options.
|
inline |
Disbale aggressive amalgamation of nodes into supernodes inside MUMPS_SYMQMAD. This is only relevant when MUMPS_SYMQAMD is enabled.
|
inline |
Do not log the assembly tree to a file. Logging of the tree is currently not supported.
|
inline |
Disable GPU off-loading.
|
inline |
Do not use the routine METIS_NodeND, but instead use the undocumented routine 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.
|
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.
|
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.
|
inline |
Disable replacement of tiny pivots.
|
inline |
When using MUMPS_SYMQAMD, enable aggressive amalgamation of nodes into supernodes. This is only used when MUMPS_SYMQAMD is enabled.
|
inline |
Log the assembly tree to a file. Currently not supported.
|
inline |
Enable off-loading to the GPU. This only works when STRUMPACK was configured with GPU support (through CUDA, MAGMA or SLATE)
|
inline |
Use the routine METIS_NodeND instead of the undocumented routine 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.
|
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.
|
inline |
Enable OpenMP tasking traversal of the supernodal tree in the sparse solver. This requires more (peak) memory, but scales better with OpenMP threads.
|
inline |
Enable positive_definite solver.
|
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).
|
inline |
Enable symmetric solver.
|
inline |
Get the GMRES restart length.
|
inline |
Returns the number of GPU streams to use.
|
inline |
Get the Gram-Schmidth orthogonalization type used in GMRES.
|
inline |
Get a reference to an object holding various options pertaining to the HODLR code, and data structures.
|
inline |
Get a (const) reference to an object holding various options pertaining to the HODLR code, and data structures.
|
inline |
Get a reference to an object holding various options pertaining to the HSS code, and data structures.
|
inline |
Get a (const) reference to an object holding various options pertaining to the HSS code, and data structures.
|
inline |
Is indirect sampling for HSS construction enabled?
|
inline |
Get the type of iterative solver to be used as outer solver.
|
inline |
Should we log the assembly tree? Currently not supported.
|
inline |
Returns the precision for lossy compression.
|
inline |
Returns the precision for lossy compression.
|
inline |
Get the matching job to use for numerical stability reordering.
|
inline |
Get the maximum number of allowed iterative solver iterations.
|
inline |
Return the nested-dissection recursion ending parameter.
|
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.
|
inline |
Get the specified nx mesh dimension.
|
inline |
Get the specified ny mesh dimension.
|
inline |
Get the specified nz mesh dimension.
|
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.
|
inline |
Info about the stats of the root front will be printed to std::cout
|
inline |
Get the type of proportional mapping to be used.
|
inline |
Get the relative tolerance to be used in the iterative solver.
|
inline |
Get the currently set fill reducing reordering method.
|
inline |
Check whether replacement of tiny pivots is enabled.
|
inline |
Get the current value of the number of additional links to include in the separator before reordering.
|
inline |
Get the currently specified width of a separator.
|
inline |
Set the absolute tolerance to use in the iterative solvers.
abs_tol | absolute tolerance |
|
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.
components | number of components per gridpoint |
|
inline |
Set the type of rank-structured compression to use.
c | compression type |
|
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().
rtol | the relative low-rank compression tolerance |
|
inline |
Set the leaf size used by any of the rank-structured formats.
|
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.
|
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.
|
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().
rtol | the relative low-rank compression tolerance |
|
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.
nx | dimension of the grid along the first axis |
ny | dimension of the grid along the second axis, should not be specified for a 1d mesh |
nz | dimension of the grid along the third axis, should not be specified for a 2d mesh |
|
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.
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.
argc | number of arguments in the argv array |
argv | list of options |
|
inline |
Set the GMRES restart length
m | GMRES restart, should be > 0 and not too crazy big |
|
inline |
Set the number of (CUDA) streams to be used in the code.
|
inline |
Set the type of Gram-Schmidt orthogonalization to use in GMRES
t | Gram-Schmidt type to use in GMRES |
|
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.
s | outer, iterative solver to use |
|
inline |
Set the accuracy for lossy compression.
|
inline |
Set the precision for lossy compression. Preferred mode is accuracy. To use precision mode, set the accuracy to a negative value. For lossless compression, set both accuracy and precision to a negative value.
|
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.
job | type of matching to perform |
|
inline |
Set the maximum number of iterations to use in any of the iterative solvers.
|
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.
|
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.
|
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.
|
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.
|
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.
|
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.
|
inline |
Print statistics, about ranks, memory etc, for the root front only.
|
inline |
Set the type of proportional mapping.
|
inline |
Set the relative tolerance to use in the iterative solvers. This will typically be the relative residual decrease.
rel_tol | relative tolerance |
|
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).
m | fill reducing reordering |
|
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.
|
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.
width | width of the separator |
|
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.
|
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.
|
inline |
Is aggressive amalgamation enabled? (only used when MUMPS_SYMQAMD is enabled)
|
inline |
Check wheter or not to use GPU off-loading.
|
inline |
Check wheter or not to use GPU-aware MPI
|
inline |
Is use of METIS_NodeND enabled? (instead of METIS_NodeNDP)
|
inline |
Is use of METIS_NodeNDP enabled? (instead of METIS_NodeND)
|
inline |
Is MUMPS_SYMQAMD enabled?
|
inline |
Check wheter or not to use OpenMP tree traversal is the sparse solver.
|
inline |
Check wheter or not to use symmetric solver.
|
inline |
Check wheter or not to use symmetric solver.
|
inline |
Check if verbose output is enabled.
|
inline |
The root front will be written to a file.