Go to the documentation of this file.
37 #include "dense/BLASLAPACKWrapper.hpp"
39 #include "BLR/BLROptions.hpp"
107 enum class EquilibrationType : char
108 { NONE=
'N', ROW=
'R', COLUMN=
'C', BOTH=
'B' };
162 {
return real_t(1.e-6); }
170 {
return real_t(1.e-10); }
174 inline int default_gpu_streams() {
return 4; }
198 using real_t =
typename RealType<scalar_t>::value_type;
217 SPOptions(
int argc,
const char*
const argv[]) : argc_(argc), argv_(argv) {
218 hss_opts_.set_verbose(
false);
219 blr_opts_.set_verbose(
false);
220 hodlr_opts_.set_verbose(
false);
255 assert(atol >=
real_t(0.));
330 assert(
nx>=1 &&
ny>=1 &&
nz>=1);
331 nx_ =
nx; ny_ =
ny; nz_ =
nz;
400 { assert(width>=1); separator_width_ = width; }
528 hss_opts_.set_rel_tol(rtol);
529 blr_opts_.set_rel_tol(rtol);
530 hodlr_opts_.set_rel_tol(rtol);
547 hss_opts_.set_abs_tol(atol);
548 blr_opts_.set_abs_tol(atol);
549 hodlr_opts_.set_abs_tol(atol);
564 hss_min_sep_size_ = s;
565 blr_min_sep_size_ = s;
566 hodlr_min_sep_size_ = s;
567 lossy_min_sep_size_ = s;
581 hss_min_front_size_ = s;
582 blr_min_front_size_ = s;
583 hodlr_min_front_size_ = s;
584 lossy_min_front_size_ = s;
593 hss_opts_.set_leaf_size(s);
594 blr_opts_.set_leaf_size(s);
595 hodlr_opts_.set_leaf_size(s);
614 { assert(l >= 0); sep_order_level_ = l; }
616 #ifndef DOXYGEN_SHOULD_SKIP_THIS
617 void enable_indirect_sampling() { indirect_sampling_ =
true; }
618 void disable_indirect_sampling() { indirect_sampling_ =
false; }
619 #endif // DOXYGEN_SHOULD_SKIP_THIS
691 int maxit()
const {
return maxit_; }
739 int nx()
const {
return nx_; }
745 int ny()
const {
return ny_; }
751 int nz()
const {
return nz_; }
820 return hss_opts_.rel_tol();
822 return blr_opts_.rel_tol();
824 return hodlr_opts_.rel_tol();
843 return hss_opts_.abs_tol();
845 return blr_opts_.abs_tol();
847 return hodlr_opts_.abs_tol();
865 return hss_min_sep_size_;
867 return blr_min_sep_size_;
869 return hodlr_min_sep_size_;
872 return lossy_min_sep_size_;
875 return std::numeric_limits<int>::max();
889 return hss_min_front_size_;
891 return blr_min_front_size_;
893 return hodlr_min_front_size_;
896 return lossy_min_front_size_;
899 return std::numeric_limits<int>::max();
914 return hss_opts_.leaf_size();
916 return blr_opts_.leaf_size();
918 return hodlr_opts_.leaf_size();
924 return std::numeric_limits<int>::max();
987 const BLR::BLROptions<scalar_t>&
BLR_options()
const {
return blr_opts_; }
1035 bool verbose_ =
true;
1038 real_t rel_tol_ = default_rel_tol<real_t>();
1039 real_t abs_tol_ = default_abs_tol<real_t>();
1041 int gmres_restart_ = 30;
1049 int components_ = 1;
1050 int separator_width_ = 1;
1051 bool use_METIS_NodeNDP_ =
false;
1052 bool use_MUMPS_SYMQAMD_ =
false;
1053 bool use_agg_amalg_ =
false;
1055 bool log_assembly_tree_ =
false;
1056 bool replace_tiny_pivots_ =
false;
1057 bool write_root_front_ =
false;
1058 bool print_root_front_stats_ =
false;
1061 bool use_gpu_ =
true;
1062 int gpu_streams_ = default_gpu_streams();
1068 int hss_min_front_size_ = 5000;
1069 int hss_min_sep_size_ = 1000;
1070 int sep_order_level_ = 1;
1071 bool indirect_sampling_ =
false;
1075 BLR::BLROptions<scalar_t> blr_opts_;
1076 int blr_min_front_size_ = 1000;
1077 int blr_min_sep_size_ = 256;
1081 int hodlr_min_front_size_ = 10000;
1082 int hodlr_min_sep_size_ = 5000;
1085 int lossy_min_front_size_ = 16;
1086 int lossy_min_sep_size_ = 8;
1087 int lossy_precision_ = 16;
1090 const char*
const* argv_ =
nullptr;
1095 #endif // SPOPTIONS_HPP
void enable_MUMPS_SYMQAMD()
Definition: StrumpackOptions.hpp:450
void set_nz(int nz)
Definition: StrumpackOptions.hpp:367
void set_nd_param(int nd_param)
Definition: StrumpackOptions.hpp:307
MatchingJob get_matching(int job)
void set_separator_width(int width)
Definition: StrumpackOptions.hpp:399
real_t abs_tol() const
Definition: StrumpackOptions.hpp:703
void enable_assembly_tree_log()
Definition: StrumpackOptions.hpp:496
int nz() const
Definition: StrumpackOptions.hpp:751
real_t rel_tol() const
Definition: StrumpackOptions.hpp:697
void disable_replace_tiny_pivots()
Definition: StrumpackOptions.hpp:645
GramSchmidtType GramSchmidt_type() const
Definition: StrumpackOptions.hpp:721
void set_components(int components)
Definition: StrumpackOptions.hpp:379
GramSchmidtType
Definition: StrumpackOptions.hpp:132
bool indirect_sampling() const
Definition: StrumpackOptions.hpp:938
int lossy_precision() const
Definition: StrumpackOptions.hpp:963
int components() const
Definition: StrumpackOptions.hpp:758
ReorderingStrategy reordering_method() const
Definition: StrumpackOptions.hpp:727
@ MAX_SMALLEST_DIAGONAL_2
int compression_leaf_size() const
Definition: StrumpackOptions.hpp:911
void disable_MUMPS_SYMQAMD()
Definition: StrumpackOptions.hpp:460
void disable_agg_amalg()
Definition: StrumpackOptions.hpp:480
KrylovSolver Krylov_solver() const
Definition: StrumpackOptions.hpp:709
HODLR::HODLROptions< scalar_t > & HODLR_options()
Definition: StrumpackOptions.hpp:1005
Class containing several options for the HSS code and data-structures.
Definition: HSSOptions.hpp:117
void disable_METIS_NodeNDP()
Definition: StrumpackOptions.hpp:422
bool log_assembly_tree() const
Definition: StrumpackOptions.hpp:801
void set_verbose(bool verbose)
Definition: StrumpackOptions.hpp:229
real_t default_rel_tol()
Definition: StrumpackOptions.hpp:161
void set_compression_min_sep_size(int s)
Definition: StrumpackOptions.hpp:562
void set_nx(int nx)
Definition: StrumpackOptions.hpp:343
BLR::BLROptions< scalar_t > & BLR_options()
Definition: StrumpackOptions.hpp:993
int compression_min_sep_size() const
Definition: StrumpackOptions.hpp:862
void set_separator_ordering_level(int l)
Definition: StrumpackOptions.hpp:613
bool use_MUMPS_SYMQAMD() const
Definition: StrumpackOptions.hpp:782
const HODLR::HODLROptions< scalar_t > & HODLR_options() const
Definition: StrumpackOptions.hpp:999
@ MAX_DIAGONAL_PRODUCT_SCALING
void set_lossy_precision(int p)
Definition: StrumpackOptions.hpp:673
const BLR::BLROptions< scalar_t > & BLR_options() const
Definition: StrumpackOptions.hpp:987
int compression_min_front_size() const
Definition: StrumpackOptions.hpp:886
Definition: StrumpackOptions.hpp:42
bool use_METIS_NodeNDP() const
Definition: StrumpackOptions.hpp:770
void set_GramSchmidt_type(GramSchmidtType t)
Definition: StrumpackOptions.hpp:285
KrylovSolver
Definition: StrumpackOptions.hpp:141
real_t default_abs_tol()
Definition: StrumpackOptions.hpp:169
void enable_METIS_NodeNDP()
Definition: StrumpackOptions.hpp:411
HSS::HSSOptions< scalar_t > & HSS_options()
Definition: StrumpackOptions.hpp:981
void set_from_command_line()
Definition: StrumpackOptions.hpp:1012
real_t compression_rel_tol() const
Definition: StrumpackOptions.hpp:817
Options for the sparse solver.
Definition: StrumpackOptions.hpp:192
void set_reordering_method(ReorderingStrategy m)
Definition: StrumpackOptions.hpp:299
void set_ny(int ny)
Definition: StrumpackOptions.hpp:355
bool use_agg_amalg() const
Definition: StrumpackOptions.hpp:789
void set_Krylov_solver(KrylovSolver s)
Definition: StrumpackOptions.hpp:271
void enable_agg_amalg()
Definition: StrumpackOptions.hpp:470
typename RealType< scalar_t >::value_type real_t
Definition: StrumpackOptions.hpp:198
void set_dimensions(int nx, int ny=1, int nz=1)
Definition: StrumpackOptions.hpp:329
void set_maxit(int maxit)
Definition: StrumpackOptions.hpp:235
SPOptions()
Definition: StrumpackOptions.hpp:204
void set_rel_tol(real_t rtol)
Definition: StrumpackOptions.hpp:243
void set_print_root_front_stats(bool b)
Definition: StrumpackOptions.hpp:679
void describe_options() const
int nd_param() const
Definition: StrumpackOptions.hpp:733
int nx() const
Definition: StrumpackOptions.hpp:739
void set_compression_rel_tol(real_t rtol)
Definition: StrumpackOptions.hpp:527
Class containing several options for the HODLR code and data-structures.
Definition: HODLROptions.hpp:115
std::string get_name(ReorderingStrategy method)
void set_gmres_restart(int m)
Definition: StrumpackOptions.hpp:278
std::string get_description(MatchingJob job)
bool use_METIS_NodeND() const
Definition: StrumpackOptions.hpp:776
real_t compression_abs_tol() const
Definition: StrumpackOptions.hpp:840
int maxit() const
Definition: StrumpackOptions.hpp:691
bool use_gpu() const
Definition: StrumpackOptions.hpp:953
void set_compression_leaf_size(int s)
Definition: StrumpackOptions.hpp:592
void set_matching(MatchingJob job)
Definition: StrumpackOptions.hpp:491
void enable_replace_tiny_pivots()
Definition: StrumpackOptions.hpp:638
void set_compression_min_front_size(int s)
Definition: StrumpackOptions.hpp:579
CompressionType compression() const
Definition: StrumpackOptions.hpp:806
int gpu_streams() const
Definition: StrumpackOptions.hpp:958
int separator_ordering_level() const
Definition: StrumpackOptions.hpp:933
Contains the HSSOptions class as well as general routines for HSS options.
void set_write_root_front(bool b)
Definition: StrumpackOptions.hpp:652
const HSS::HSSOptions< scalar_t > & HSS_options() const
Definition: StrumpackOptions.hpp:975
bool verbose() const
Definition: StrumpackOptions.hpp:685
void disable_METIS_NodeND()
Definition: StrumpackOptions.hpp:441
bool replace_tiny_pivots() const
Definition: StrumpackOptions.hpp:943
void disable_gpu()
Definition: StrumpackOptions.hpp:663
void set_gpu_streams(int s)
Definition: StrumpackOptions.hpp:668
void enable_METIS_NodeND()
Definition: StrumpackOptions.hpp:432
bool write_root_front() const
Definition: StrumpackOptions.hpp:948
int separator_width() const
Definition: StrumpackOptions.hpp:764
bool is_parallel(ReorderingStrategy method)
void set_compression_abs_tol(real_t atol)
Definition: StrumpackOptions.hpp:546
Contains the class holding HODLR matrix options.
int gmres_restart() const
Definition: StrumpackOptions.hpp:715
void enable_gpu()
Definition: StrumpackOptions.hpp:658
int ny() const
Definition: StrumpackOptions.hpp:745
SPOptions(int argc, const char *const argv[])
Definition: StrumpackOptions.hpp:217
bool print_root_front_stats() const
Definition: StrumpackOptions.hpp:969
MatchingJob matching() const
Definition: StrumpackOptions.hpp:795
void disable_assembly_tree_log()
Definition: StrumpackOptions.hpp:502
CompressionType
Definition: StrumpackOptions.hpp:74
void set_abs_tol(real_t atol)
Definition: StrumpackOptions.hpp:254
ReorderingStrategy
Definition: StrumpackOptions.hpp:48
void set_compression(CompressionType c)
Definition: StrumpackOptions.hpp:512
MatchingJob
Definition: StrumpackOptions.hpp:95