Go to the documentation of this file.
37 #include "dense/BLASLAPACKWrapper.hpp"
39 #include "BLR/BLROptions.hpp"
110 enum class EquilibrationType : char
111 { NONE=
'N', ROW=
'R', COLUMN=
'C', BOTH=
'B' };
165 {
return real_t(1.e-6); }
173 {
return real_t(1.e-10); }
177 inline int default_gpu_streams() {
return 4; }
201 using real_t =
typename RealType<scalar_t>::value_type;
220 SPOptions(
int argc,
const char*
const argv[]) : argc_(argc), argv_(argv) {
221 hss_opts_.set_verbose(
false);
222 blr_opts_.set_verbose(
false);
223 hodlr_opts_.set_verbose(
false);
258 assert(atol >=
real_t(0.));
342 assert(
nx>=1 &&
ny>=1 &&
nz>=1);
343 nx_ =
nx; ny_ =
ny; nz_ =
nz;
412 { assert(width>=1); separator_width_ = width; }
540 hss_opts_.set_rel_tol(rtol);
541 blr_opts_.set_rel_tol(rtol);
542 hodlr_opts_.set_rel_tol(rtol);
559 hss_opts_.set_abs_tol(atol);
560 blr_opts_.set_abs_tol(atol);
561 hodlr_opts_.set_abs_tol(atol);
576 hss_min_sep_size_ = s;
577 blr_min_sep_size_ = s;
578 hodlr_min_sep_size_ = s;
579 lossy_min_sep_size_ = s;
581 void set_hss_min_sep_size(
int s) {
583 hss_min_sep_size_ = s;
585 void set_hodlr_min_sep_size(
int s) {
587 hodlr_min_sep_size_ = s;
589 void set_blr_min_sep_size(
int s) {
591 blr_min_sep_size_ = s;
593 void set_lossy_min_sep_size(
int s) {
595 lossy_min_sep_size_ = s;
609 hss_min_front_size_ = s;
610 blr_min_front_size_ = s;
611 hodlr_min_front_size_ = s;
612 lossy_min_front_size_ = s;
614 void set_hss_min_front_size(
int s) {
616 hss_min_front_size_ = s;
618 void set_hodlr_min_front_size(
int s) {
620 hodlr_min_front_size_ = s;
622 void set_blr_min_front_size(
int s) {
624 blr_min_front_size_ = s;
626 void set_lossy_min_front_size(
int s) {
628 lossy_min_front_size_ = s;
637 hss_opts_.set_leaf_size(s);
638 blr_opts_.set_leaf_size(s);
639 hodlr_opts_.set_leaf_size(s);
658 { assert(l >= 0); sep_order_level_ = l; }
660 #ifndef DOXYGEN_SHOULD_SKIP_THIS
661 void enable_indirect_sampling() { indirect_sampling_ =
true; }
662 void disable_indirect_sampling() { indirect_sampling_ =
false; }
663 #endif // DOXYGEN_SHOULD_SKIP_THIS
746 int maxit()
const {
return maxit_; }
802 int nx()
const {
return nx_; }
808 int ny()
const {
return ny_; }
814 int nz()
const {
return nz_; }
883 return hss_opts_.rel_tol();
885 return blr_opts_.rel_tol();
887 return hodlr_opts_.rel_tol();
889 if (l==0)
return hodlr_opts_.rel_tol();
890 else return blr_opts_.rel_tol();
909 return hss_opts_.abs_tol();
911 return blr_opts_.abs_tol();
913 return hodlr_opts_.abs_tol();
915 if (l==0)
return hodlr_opts_.abs_tol();
916 else return blr_opts_.abs_tol();
934 return hss_min_sep_size_;
936 return blr_min_sep_size_;
938 return hodlr_min_sep_size_;
940 if (l==0)
return hodlr_min_sep_size_;
941 else return blr_min_sep_size_;
944 return lossy_min_sep_size_;
947 return std::numeric_limits<int>::max();
950 int hss_min_sep_size()
const {
951 return hss_min_sep_size_;
953 int hodlr_min_sep_size()
const {
954 return hodlr_min_sep_size_;
956 int blr_min_sep_size()
const {
957 return blr_min_sep_size_;
959 int lossy_min_sep_size()
const {
960 return lossy_min_sep_size_;
973 return hss_min_front_size_;
975 return blr_min_front_size_;
977 return hodlr_min_front_size_;
979 if (l==0)
return hodlr_min_front_size_;
980 else return blr_min_front_size_;
983 return lossy_min_front_size_;
986 return std::numeric_limits<int>::max();
989 int hss_min_front_size()
const {
990 return hss_min_front_size_;
992 int hodlr_min_front_size()
const {
993 return hodlr_min_front_size_;
995 int blr_min_front_size()
const {
996 return blr_min_front_size_;
998 int lossy_min_front_size()
const {
999 return lossy_min_front_size_;
1013 return hss_opts_.leaf_size();
1015 return blr_opts_.leaf_size();
1017 return hodlr_opts_.leaf_size();
1019 if (l==0)
return hodlr_opts_.leaf_size();
1020 else return blr_opts_.leaf_size();
1026 return std::numeric_limits<int>::max();
1101 const BLR::BLROptions<scalar_t>&
BLR_options()
const {
return blr_opts_; }
1149 bool verbose_ =
true;
1152 real_t rel_tol_ = default_rel_tol<real_t>();
1153 real_t abs_tol_ = default_abs_tol<real_t>();
1155 int gmres_restart_ = 30;
1159 int nd_planar_levels_ = 0;
1164 int components_ = 1;
1165 int separator_width_ = 1;
1166 bool use_METIS_NodeNDP_ =
false;
1167 bool use_MUMPS_SYMQAMD_ =
false;
1168 bool use_agg_amalg_ =
false;
1170 bool log_assembly_tree_ =
false;
1171 bool replace_tiny_pivots_ =
false;
1172 real_t pivot_ = std::sqrt(blas::lamch<real_t>(
'E'));
1173 bool write_root_front_ =
false;
1174 bool print_root_front_stats_ =
false;
1177 bool use_gpu_ =
true;
1178 int gpu_streams_ = default_gpu_streams();
1184 int hss_min_front_size_ = 5000;
1185 int hss_min_sep_size_ = 1000;
1186 int sep_order_level_ = 1;
1187 bool indirect_sampling_ =
false;
1191 BLR::BLROptions<scalar_t> blr_opts_;
1192 int blr_min_front_size_ = 1000;
1193 int blr_min_sep_size_ = 256;
1197 int hodlr_min_front_size_ = 10000;
1198 int hodlr_min_sep_size_ = 5000;
1201 int lossy_min_front_size_ = 16;
1202 int lossy_min_sep_size_ = 8;
1203 int lossy_precision_ = 16;
1206 const char*
const* argv_ =
nullptr;
1211 #endif // SPOPTIONS_HPP
void enable_MUMPS_SYMQAMD()
Definition: StrumpackOptions.hpp:462
void set_nz(int nz)
Definition: StrumpackOptions.hpp:379
void set_nd_param(int nd_param)
Definition: StrumpackOptions.hpp:310
MatchingJob get_matching(int job)
void set_separator_width(int width)
Definition: StrumpackOptions.hpp:411
real_t abs_tol() const
Definition: StrumpackOptions.hpp:758
void enable_assembly_tree_log()
Definition: StrumpackOptions.hpp:508
int nz() const
Definition: StrumpackOptions.hpp:814
real_t rel_tol() const
Definition: StrumpackOptions.hpp:752
void disable_replace_tiny_pivots()
Definition: StrumpackOptions.hpp:689
void set_nd_planar_levels(int nd_planar_levels)
Definition: StrumpackOptions.hpp:319
GramSchmidtType GramSchmidt_type() const
Definition: StrumpackOptions.hpp:776
void set_components(int components)
Definition: StrumpackOptions.hpp:391
GramSchmidtType
Definition: StrumpackOptions.hpp:135
bool indirect_sampling() const
Definition: StrumpackOptions.hpp:1040
int lossy_precision() const
Definition: StrumpackOptions.hpp:1077
int components() const
Definition: StrumpackOptions.hpp:821
ReorderingStrategy reordering_method() const
Definition: StrumpackOptions.hpp:782
@ MAX_SMALLEST_DIAGONAL_2
void disable_MUMPS_SYMQAMD()
Definition: StrumpackOptions.hpp:472
void disable_agg_amalg()
Definition: StrumpackOptions.hpp:492
KrylovSolver Krylov_solver() const
Definition: StrumpackOptions.hpp:764
HODLR::HODLROptions< scalar_t > & HODLR_options()
Definition: StrumpackOptions.hpp:1119
Class containing several options for the HSS code and data-structures.
Definition: HSSOptions.hpp:117
void disable_METIS_NodeNDP()
Definition: StrumpackOptions.hpp:434
bool log_assembly_tree() const
Definition: StrumpackOptions.hpp:864
void set_verbose(bool verbose)
Definition: StrumpackOptions.hpp:232
real_t default_rel_tol()
Definition: StrumpackOptions.hpp:164
void set_compression_min_sep_size(int s)
Definition: StrumpackOptions.hpp:574
void set_nx(int nx)
Definition: StrumpackOptions.hpp:355
BLR::BLROptions< scalar_t > & BLR_options()
Definition: StrumpackOptions.hpp:1107
real_t compression_rel_tol(int l=0) const
Definition: StrumpackOptions.hpp:880
void set_separator_ordering_level(int l)
Definition: StrumpackOptions.hpp:657
bool use_MUMPS_SYMQAMD() const
Definition: StrumpackOptions.hpp:845
const HODLR::HODLROptions< scalar_t > & HODLR_options() const
Definition: StrumpackOptions.hpp:1113
@ MAX_DIAGONAL_PRODUCT_SCALING
void set_lossy_precision(int p)
Definition: StrumpackOptions.hpp:728
const BLR::BLROptions< scalar_t > & BLR_options() const
Definition: StrumpackOptions.hpp:1101
Definition: StrumpackOptions.hpp:42
bool use_METIS_NodeNDP() const
Definition: StrumpackOptions.hpp:833
void set_GramSchmidt_type(GramSchmidtType t)
Definition: StrumpackOptions.hpp:288
KrylovSolver
Definition: StrumpackOptions.hpp:144
real_t default_abs_tol()
Definition: StrumpackOptions.hpp:172
void enable_METIS_NodeNDP()
Definition: StrumpackOptions.hpp:423
HSS::HSSOptions< scalar_t > & HSS_options()
Definition: StrumpackOptions.hpp:1095
void set_from_command_line()
Definition: StrumpackOptions.hpp:1126
Options for the sparse solver.
Definition: StrumpackOptions.hpp:195
void set_reordering_method(ReorderingStrategy m)
Definition: StrumpackOptions.hpp:302
void set_ny(int ny)
Definition: StrumpackOptions.hpp:367
bool use_agg_amalg() const
Definition: StrumpackOptions.hpp:852
void set_Krylov_solver(KrylovSolver s)
Definition: StrumpackOptions.hpp:274
void enable_agg_amalg()
Definition: StrumpackOptions.hpp:482
typename RealType< refine_t >::value_type real_t
Definition: StrumpackOptions.hpp:201
int compression_min_front_size(int l=0) const
Definition: StrumpackOptions.hpp:970
void set_dimensions(int nx, int ny=1, int nz=1)
Definition: StrumpackOptions.hpp:341
void set_maxit(int maxit)
Definition: StrumpackOptions.hpp:238
SPOptions()
Definition: StrumpackOptions.hpp:207
void set_rel_tol(real_t rtol)
Definition: StrumpackOptions.hpp:246
void set_print_root_front_stats(bool b)
Definition: StrumpackOptions.hpp:734
void describe_options() const
int nd_param() const
Definition: StrumpackOptions.hpp:788
int nx() const
Definition: StrumpackOptions.hpp:802
void set_compression_rel_tol(real_t rtol)
Definition: StrumpackOptions.hpp:539
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:281
std::string get_description(MatchingJob job)
bool use_METIS_NodeND() const
Definition: StrumpackOptions.hpp:839
real_t compression_abs_tol(int l=0) const
Definition: StrumpackOptions.hpp:906
int maxit() const
Definition: StrumpackOptions.hpp:746
bool use_gpu() const
Definition: StrumpackOptions.hpp:1067
void set_compression_leaf_size(int s)
Definition: StrumpackOptions.hpp:636
void set_matching(MatchingJob job)
Definition: StrumpackOptions.hpp:503
void enable_replace_tiny_pivots()
Definition: StrumpackOptions.hpp:682
void set_compression_min_front_size(int s)
Definition: StrumpackOptions.hpp:607
CompressionType compression() const
Definition: StrumpackOptions.hpp:869
int gpu_streams() const
Definition: StrumpackOptions.hpp:1072
int separator_ordering_level() const
Definition: StrumpackOptions.hpp:1035
Contains the HSSOptions class as well as general routines for HSS options.
void set_write_root_front(bool b)
Definition: StrumpackOptions.hpp:707
int compression_leaf_size(int l=0) const
Definition: StrumpackOptions.hpp:1010
const HSS::HSSOptions< scalar_t > & HSS_options() const
Definition: StrumpackOptions.hpp:1089
bool verbose() const
Definition: StrumpackOptions.hpp:740
void disable_METIS_NodeND()
Definition: StrumpackOptions.hpp:453
bool replace_tiny_pivots() const
Definition: StrumpackOptions.hpp:1047
void disable_gpu()
Definition: StrumpackOptions.hpp:718
void set_gpu_streams(int s)
Definition: StrumpackOptions.hpp:723
void enable_METIS_NodeND()
Definition: StrumpackOptions.hpp:444
void set_pivot_threshold(real_t thresh)
Definition: StrumpackOptions.hpp:700
bool write_root_front() const
Definition: StrumpackOptions.hpp:1062
int separator_width() const
Definition: StrumpackOptions.hpp:827
real_t pivot_threshold() const
Definition: StrumpackOptions.hpp:1057
bool is_parallel(ReorderingStrategy method)
void set_compression_abs_tol(real_t atol)
Definition: StrumpackOptions.hpp:558
Contains the class holding HODLR matrix options.
int gmres_restart() const
Definition: StrumpackOptions.hpp:770
void enable_gpu()
Definition: StrumpackOptions.hpp:713
int ny() const
Definition: StrumpackOptions.hpp:808
int nd_planar_levels() const
Definition: StrumpackOptions.hpp:796
SPOptions(int argc, const char *const argv[])
Definition: StrumpackOptions.hpp:220
bool print_root_front_stats() const
Definition: StrumpackOptions.hpp:1083
MatchingJob matching() const
Definition: StrumpackOptions.hpp:858
void disable_assembly_tree_log()
Definition: StrumpackOptions.hpp:514
CompressionType
Definition: StrumpackOptions.hpp:74
void set_abs_tol(real_t atol)
Definition: StrumpackOptions.hpp:257
int compression_min_sep_size(int l=0) const
Definition: StrumpackOptions.hpp:931
ReorderingStrategy
Definition: StrumpackOptions.hpp:48
void set_compression(CompressionType c)
Definition: StrumpackOptions.hpp:524
MatchingJob
Definition: StrumpackOptions.hpp:98