StrumpackOptions.hpp
Go to the documentation of this file.
1 /*
2  * STRUMPACK -- STRUctured Matrices PACKage, Copyright (c) 2014, The
3  * Regents of the University of California, through Lawrence Berkeley
4  * National Laboratory (subject to receipt of any required approvals
5  * from the U.S. Dept. of Energy). All rights reserved.
6  *
7  * If you have questions about your rights to use or distribute this
8  * software, please contact Berkeley Lab's Technology Transfer
9  * Department at TTD@lbl.gov.
10  *
11  * NOTICE. This software is owned by the U.S. Department of Energy. As
12  * such, the U.S. Government has been granted for itself and others
13  * acting on its behalf a paid-up, nonexclusive, irrevocable,
14  * worldwide license in the Software to reproduce, prepare derivative
15  * works, and perform publicly and display publicly. Beginning five
16  * (5) years after the date permission to assert copyright is obtained
17  * from the U.S. Department of Energy, and subject to any subsequent
18  * five (5) year renewals, the U.S. Government is granted for itself
19  * and others acting on its behalf a paid-up, nonexclusive,
20  * irrevocable, worldwide license in the Software to reproduce,
21  * prepare derivative works, distribute copies to the public, perform
22  * publicly and display publicly, and to permit others to do so.
23  *
24  * Developers: Pieter Ghysels, Francois-Henry Rouet, Xiaoye S. Li.
25  * (Lawrence Berkeley National Lab, Computational Research
26  * Division).
27  */
32 #ifndef SPOPTIONS_HPP
33 #define SPOPTIONS_HPP
34 
35 #include <limits>
36 
37 #include "dense/BLASLAPACKWrapper.hpp"
38 #include "HSS/HSSOptions.hpp"
39 #include "BLR/BLROptions.hpp"
40 #include "HODLR/HODLROptions.hpp"
41 
42 namespace strumpack {
43 
48  enum class ReorderingStrategy {
49  NATURAL,
50  METIS,
51  PARMETIS,
52  SCOTCH,
53  PTSCOTCH,
54  RCM,
55  GEOMETRIC
57  };
58 
62  std::string get_name(ReorderingStrategy method);
63 
67  bool is_parallel(ReorderingStrategy method);
68 
74  enum class CompressionType {
75  NONE,
76  HSS,
77  BLR,
78  HODLR,
80  BLR_HODLR,
83  LOSSLESS,
84  LOSSY
85  };
86 
90  std::string get_name(CompressionType comp);
91 
92 
98  enum class MatchingJob {
99  NONE,
107  COMBBLAS
108  };
109 
110  enum class EquilibrationType : char
111  { NONE='N', ROW='R', COLUMN='C', BOTH='B' };
112 
113 
117  MatchingJob get_matching(int job);
118 
123  int get_matching(MatchingJob job);
124 
128  std::string get_description(MatchingJob job);
129 
130 
135  enum class GramSchmidtType {
136  CLASSICAL,
137  MODIFIED
138  };
139 
144  enum class KrylovSolver {
145  AUTO,
147  DIRECT,
149  REFINE,
150  PREC_GMRES,
152  GMRES,
153  PREC_BICGSTAB,
155  BICGSTAB
156  };
157 
164  template<typename real_t> inline real_t default_rel_tol()
165  { return real_t(1.e-6); }
172  template<typename real_t> inline real_t default_abs_tol()
173  { return real_t(1.e-10); }
174  template<> inline float default_rel_tol() { return 1.e-4; }
175  template<> inline float default_abs_tol() { return 1.e-6; }
176 
177  inline int default_gpu_streams() { return 4; }
178 
195  template<typename scalar_t> class SPOptions {
196  public:
201  using real_t = typename RealType<scalar_t>::value_type;
202 
207  SPOptions() { hss_opts_.set_verbose(false); }
208 
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);
224  }
225 
232  void set_verbose(bool verbose) { verbose_ = verbose; }
233 
238  void set_maxit(int maxit) { assert(maxit >= 1); maxit_ = maxit; }
239 
246  void set_rel_tol(real_t rtol) {
247  assert(rtol <= real_t(1.) && rtol >= real_t(0.));
248  rel_tol_ = rtol;
249  }
250 
257  void set_abs_tol(real_t atol) {
258  assert(atol >= real_t(0.));
259  abs_tol_ = atol;
260  }
261 
274  void set_Krylov_solver(KrylovSolver s) { Krylov_solver_ = s; }
275 
281  void set_gmres_restart(int m) { assert(m >= 1); gmres_restart_ = m; }
282 
288  void set_GramSchmidt_type(GramSchmidtType t) { Gram_Schmidt_type_ = t; }
289 
302  void set_reordering_method(ReorderingStrategy m) { reordering_method_ = m; }
303 
311  { assert(nd_param>=0); nd_param_ = nd_param; }
312 
320  { assert(nd_planar_levels>=0); nd_planar_levels_ = nd_planar_levels; }
321 
341  void set_dimensions(int nx, int ny=1, int nz=1) {
342  assert(nx>=1 && ny>=1 && nz>=1);
343  nx_ = nx; ny_ = ny; nz_ = nz;
344  }
345 
355  void set_nx(int nx) {assert(nx>=1); nx_ = nx; }
356 
367  void set_ny(int ny) {assert(ny>=1); ny_ = ny; }
368 
379  void set_nz(int nz) { assert(nz>=1); nz_ = nz; }
380 
392  { assert(components>=1); components_ = components; }
393 
411  void set_separator_width(int width)
412  { assert(width>=1); separator_width_ = width; }
413 
423  void enable_METIS_NodeNDP() { use_METIS_NodeNDP_ = true; }
424 
434  void disable_METIS_NodeNDP() { use_METIS_NodeNDP_ = false; }
435 
436 
444  void enable_METIS_NodeND() { use_METIS_NodeNDP_ = false; }
445 
453  void disable_METIS_NodeND() { use_METIS_NodeNDP_ = true; }
454 
462  void enable_MUMPS_SYMQAMD() { use_MUMPS_SYMQAMD_ = true; }
463 
472  void disable_MUMPS_SYMQAMD() { use_MUMPS_SYMQAMD_ = false; }
473 
482  void enable_agg_amalg() { use_agg_amalg_ = true; }
483 
492  void disable_agg_amalg() { use_agg_amalg_ = false; }
493 
503  void set_matching(MatchingJob job) { matching_job_ = job; }
504 
508  void enable_assembly_tree_log() { log_assembly_tree_ = true; }
509 
514  void disable_assembly_tree_log() { log_assembly_tree_ = false; }
515 
524  void set_compression(CompressionType c) { comp_ = c; }
525 
540  hss_opts_.set_rel_tol(rtol);
541  blr_opts_.set_rel_tol(rtol);
542  hodlr_opts_.set_rel_tol(rtol);
543  }
544 
559  hss_opts_.set_abs_tol(atol);
560  blr_opts_.set_abs_tol(atol);
561  hodlr_opts_.set_abs_tol(atol);
562  }
563 
575  assert(s >= 0);
576  hss_min_sep_size_ = s;
577  blr_min_sep_size_ = s;
578  hodlr_min_sep_size_ = s;
579  lossy_min_sep_size_ = s;
580  }
581  void set_hss_min_sep_size(int s) {
582  assert(s >= 0);
583  hss_min_sep_size_ = s;
584  }
585  void set_hodlr_min_sep_size(int s) {
586  assert(s >= 0);
587  hodlr_min_sep_size_ = s;
588  }
589  void set_blr_min_sep_size(int s) {
590  assert(s >= 0);
591  blr_min_sep_size_ = s;
592  }
593  void set_lossy_min_sep_size(int s) {
594  assert(s >= 0);
595  lossy_min_sep_size_ = s;
596  }
597 
608  assert(s >= 0);
609  hss_min_front_size_ = s;
610  blr_min_front_size_ = s;
611  hodlr_min_front_size_ = s;
612  lossy_min_front_size_ = s;
613  }
614  void set_hss_min_front_size(int s) {
615  assert(s >= 0);
616  hss_min_front_size_ = s;
617  }
618  void set_hodlr_min_front_size(int s) {
619  assert(s >= 0);
620  hodlr_min_front_size_ = s;
621  }
622  void set_blr_min_front_size(int s) {
623  assert(s >= 0);
624  blr_min_front_size_ = s;
625  }
626  void set_lossy_min_front_size(int s) {
627  assert(s >= 0);
628  lossy_min_front_size_ = s;
629  }
630 
637  hss_opts_.set_leaf_size(s);
638  blr_opts_.set_leaf_size(s);
639  hodlr_opts_.set_leaf_size(s);
640  }
641 
658  { assert(l >= 0); sep_order_level_ = l; }
659 
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
664 
682  void enable_replace_tiny_pivots() { replace_tiny_pivots_ = true; }
683 
689  void disable_replace_tiny_pivots() { replace_tiny_pivots_ = false; }
690 
700  void set_pivot_threshold(real_t thresh) { pivot_ = thresh; }
701 
707  void set_write_root_front(bool b) { write_root_front_ = b; }
708 
713  void enable_gpu() { use_gpu_ = true; }
714 
718  void disable_gpu() { use_gpu_ = false; }
719 
723  void set_gpu_streams(int s) { gpu_streams_ = s; }
724 
728  void set_lossy_precision(int p) { lossy_precision_ = p; }
729 
734  void set_print_root_front_stats(bool b) { print_root_front_stats_ = b; }
735 
740  bool verbose() const { return verbose_; }
741 
746  int maxit() const { return maxit_; }
747 
752  real_t rel_tol() const { return rel_tol_; }
753 
758  real_t abs_tol() const { return abs_tol_; }
759 
764  KrylovSolver Krylov_solver() const { return Krylov_solver_; }
765 
770  int gmres_restart() const { return gmres_restart_; }
771 
776  GramSchmidtType GramSchmidt_type() const { return Gram_Schmidt_type_; }
777 
782  ReorderingStrategy reordering_method() const { return reordering_method_; }
783 
788  int nd_param() const { return nd_param_; }
789 
796  int nd_planar_levels() const { return nd_planar_levels_; }
797 
802  int nx() const { return nx_; }
803 
808  int ny() const { return ny_; }
809 
814  int nz() const { return nz_; }
815 
821  int components() const { return components_; }
822 
827  int separator_width() const { return separator_width_; }
828 
833  bool use_METIS_NodeNDP() const { return use_METIS_NodeNDP_; }
834 
839  bool use_METIS_NodeND() const { return !use_METIS_NodeNDP_; }
840 
845  bool use_MUMPS_SYMQAMD() const { return use_MUMPS_SYMQAMD_; }
846 
852  bool use_agg_amalg() const { return use_agg_amalg_; }
853 
858  MatchingJob matching() const { return matching_job_; }
859 
864  bool log_assembly_tree() const { return log_assembly_tree_; }
865 
869  CompressionType compression() const { return comp_; }
870 
871 
880  real_t compression_rel_tol(int l=0) const {
881  switch (comp_) {
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();
894  default: return 0.;
895  }
896  }
897 
906  real_t compression_abs_tol(int l=0) const {
907  switch (comp_) {
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();
920  default: return 0.;
921  }
922  }
923 
931  int compression_min_sep_size(int l=0) const {
932  switch (comp_) {
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_;
946  default:
947  return std::numeric_limits<int>::max();
948  }
949  }
950  int hss_min_sep_size() const {
951  return hss_min_sep_size_;
952  }
953  int hodlr_min_sep_size() const {
954  return hodlr_min_sep_size_;
955  }
956  int blr_min_sep_size() const {
957  return blr_min_sep_size_;
958  }
959  int lossy_min_sep_size() const {
960  return lossy_min_sep_size_;
961  }
962 
970  int compression_min_front_size(int l=0) const {
971  switch (comp_) {
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_;
985  default:
986  return std::numeric_limits<int>::max();
987  }
988  }
989  int hss_min_front_size() const {
990  return hss_min_front_size_;
991  }
992  int hodlr_min_front_size() const {
993  return hodlr_min_front_size_;
994  }
995  int blr_min_front_size() const {
996  return blr_min_front_size_;
997  }
998  int lossy_min_front_size() const {
999  return lossy_min_front_size_;
1000  }
1001 
1010  int compression_leaf_size(int l=0) const {
1011  switch (comp_) {
1012  case CompressionType::HSS:
1013  return hss_opts_.leaf_size();
1014  case CompressionType::BLR:
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();
1023  return 4;
1024  case CompressionType::NONE:
1025  default:
1026  return std::numeric_limits<int>::max();
1027  }
1028  }
1029 
1035  int separator_ordering_level() const { return sep_order_level_; }
1036 
1040  bool indirect_sampling() const { return indirect_sampling_; }
1041 
1047  bool replace_tiny_pivots() const { return replace_tiny_pivots_; }
1048 
1057  real_t pivot_threshold() const { return pivot_; }
1058 
1062  bool write_root_front() const { return write_root_front_; }
1063 
1067  bool use_gpu() const { return use_gpu_; }
1068 
1072  int gpu_streams() const { return gpu_streams_; }
1073 
1077  int lossy_precision() const { return lossy_precision_; }
1078 
1083  bool print_root_front_stats() const { return print_root_front_stats_; }
1084 
1089  const HSS::HSSOptions<scalar_t>& HSS_options() const { return hss_opts_; }
1090 
1095  HSS::HSSOptions<scalar_t>& HSS_options() { return hss_opts_; }
1096 
1101  const BLR::BLROptions<scalar_t>& BLR_options() const { return blr_opts_; }
1102 
1107  BLR::BLROptions<scalar_t>& BLR_options() { return blr_opts_; }
1108 
1113  const HODLR::HODLROptions<scalar_t>& HODLR_options() const { return hodlr_opts_; }
1114 
1120 
1127 
1140  void set_from_command_line(int argc, const char* const* cargv);
1141 
1146  void describe_options() const;
1147 
1148  private:
1149  bool verbose_ = true;
1151  int maxit_ = 5000;
1152  real_t rel_tol_ = default_rel_tol<real_t>();
1153  real_t abs_tol_ = default_abs_tol<real_t>();
1154  KrylovSolver Krylov_solver_ = KrylovSolver::AUTO;
1155  int gmres_restart_ = 30;
1156  GramSchmidtType Gram_Schmidt_type_ = GramSchmidtType::MODIFIED;
1158  ReorderingStrategy reordering_method_ = ReorderingStrategy::METIS;
1159  int nd_planar_levels_ = 0;
1160  int nd_param_ = 8;
1161  int nx_ = 1;
1162  int ny_ = 1;
1163  int nz_ = 1;
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;
1175 
1177  bool use_gpu_ = true;
1178  int gpu_streams_ = default_gpu_streams();
1179 
1182 
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;
1188  HSS::HSSOptions<scalar_t> hss_opts_;
1189 
1191  BLR::BLROptions<scalar_t> blr_opts_;
1192  int blr_min_front_size_ = 1000;
1193  int blr_min_sep_size_ = 256;
1194 
1196  HODLR::HODLROptions<scalar_t> hodlr_opts_;
1197  int hodlr_min_front_size_ = 10000;
1198  int hodlr_min_sep_size_ = 5000;
1199 
1201  int lossy_min_front_size_ = 16;
1202  int lossy_min_sep_size_ = 8;
1203  int lossy_precision_ = 16;
1204 
1205  int argc_ = 0;
1206  const char* const* argv_ = nullptr;
1207  };
1208 
1209 } // end namespace strumpack
1210 
1211 #endif // SPOPTIONS_HPP
strumpack::SPOptions::enable_MUMPS_SYMQAMD
void enable_MUMPS_SYMQAMD()
Definition: StrumpackOptions.hpp:462
strumpack::KrylovSolver::PREC_GMRES
@ PREC_GMRES
strumpack::SPOptions::set_nz
void set_nz(int nz)
Definition: StrumpackOptions.hpp:379
strumpack::ReorderingStrategy::NATURAL
@ NATURAL
strumpack::SPOptions::set_nd_param
void set_nd_param(int nd_param)
Definition: StrumpackOptions.hpp:310
strumpack::get_matching
MatchingJob get_matching(int job)
strumpack::SPOptions::set_separator_width
void set_separator_width(int width)
Definition: StrumpackOptions.hpp:411
strumpack::SPOptions::abs_tol
real_t abs_tol() const
Definition: StrumpackOptions.hpp:758
strumpack::SPOptions::enable_assembly_tree_log
void enable_assembly_tree_log()
Definition: StrumpackOptions.hpp:508
strumpack::SPOptions::nz
int nz() const
Definition: StrumpackOptions.hpp:814
strumpack::SPOptions::rel_tol
real_t rel_tol() const
Definition: StrumpackOptions.hpp:752
strumpack::SPOptions::disable_replace_tiny_pivots
void disable_replace_tiny_pivots()
Definition: StrumpackOptions.hpp:689
strumpack::SPOptions::set_nd_planar_levels
void set_nd_planar_levels(int nd_planar_levels)
Definition: StrumpackOptions.hpp:319
strumpack::SPOptions::GramSchmidt_type
GramSchmidtType GramSchmidt_type() const
Definition: StrumpackOptions.hpp:776
strumpack::SPOptions::set_components
void set_components(int components)
Definition: StrumpackOptions.hpp:391
strumpack::GramSchmidtType
GramSchmidtType
Definition: StrumpackOptions.hpp:135
strumpack::SPOptions::indirect_sampling
bool indirect_sampling() const
Definition: StrumpackOptions.hpp:1040
strumpack::KrylovSolver::GMRES
@ GMRES
strumpack::SPOptions::lossy_precision
int lossy_precision() const
Definition: StrumpackOptions.hpp:1077
strumpack::SPOptions::components
int components() const
Definition: StrumpackOptions.hpp:821
strumpack::SPOptions::reordering_method
ReorderingStrategy reordering_method() const
Definition: StrumpackOptions.hpp:782
strumpack::MatchingJob::MAX_SMALLEST_DIAGONAL_2
@ MAX_SMALLEST_DIAGONAL_2
strumpack::SPOptions::disable_MUMPS_SYMQAMD
void disable_MUMPS_SYMQAMD()
Definition: StrumpackOptions.hpp:472
strumpack::SPOptions::disable_agg_amalg
void disable_agg_amalg()
Definition: StrumpackOptions.hpp:492
strumpack::SPOptions::Krylov_solver
KrylovSolver Krylov_solver() const
Definition: StrumpackOptions.hpp:764
strumpack::CompressionType::LOSSLESS
@ LOSSLESS
strumpack::SPOptions::HODLR_options
HODLR::HODLROptions< scalar_t > & HODLR_options()
Definition: StrumpackOptions.hpp:1119
strumpack::GramSchmidtType::CLASSICAL
@ CLASSICAL
strumpack::HSS::HSSOptions
Class containing several options for the HSS code and data-structures.
Definition: HSSOptions.hpp:117
strumpack::SPOptions::disable_METIS_NodeNDP
void disable_METIS_NodeNDP()
Definition: StrumpackOptions.hpp:434
strumpack::KrylovSolver::AUTO
@ AUTO
strumpack::ReorderingStrategy::PTSCOTCH
@ PTSCOTCH
strumpack::SPOptions::log_assembly_tree
bool log_assembly_tree() const
Definition: StrumpackOptions.hpp:864
strumpack::SPOptions::set_verbose
void set_verbose(bool verbose)
Definition: StrumpackOptions.hpp:232
strumpack::default_rel_tol
real_t default_rel_tol()
Definition: StrumpackOptions.hpp:164
strumpack::MatchingJob::MAX_DIAGONAL_SUM
@ MAX_DIAGONAL_SUM
strumpack::SPOptions::set_compression_min_sep_size
void set_compression_min_sep_size(int s)
Definition: StrumpackOptions.hpp:574
strumpack::SPOptions::set_nx
void set_nx(int nx)
Definition: StrumpackOptions.hpp:355
strumpack::SPOptions::BLR_options
BLR::BLROptions< scalar_t > & BLR_options()
Definition: StrumpackOptions.hpp:1107
strumpack::SPOptions::compression_rel_tol
real_t compression_rel_tol(int l=0) const
Definition: StrumpackOptions.hpp:880
strumpack::SPOptions::set_separator_ordering_level
void set_separator_ordering_level(int l)
Definition: StrumpackOptions.hpp:657
strumpack::SPOptions::use_MUMPS_SYMQAMD
bool use_MUMPS_SYMQAMD() const
Definition: StrumpackOptions.hpp:845
strumpack::KrylovSolver::DIRECT
@ DIRECT
strumpack::SPOptions::HODLR_options
const HODLR::HODLROptions< scalar_t > & HODLR_options() const
Definition: StrumpackOptions.hpp:1113
strumpack::MatchingJob::MAX_DIAGONAL_PRODUCT_SCALING
@ MAX_DIAGONAL_PRODUCT_SCALING
strumpack::SPOptions::set_lossy_precision
void set_lossy_precision(int p)
Definition: StrumpackOptions.hpp:728
strumpack::SPOptions::BLR_options
const BLR::BLROptions< scalar_t > & BLR_options() const
Definition: StrumpackOptions.hpp:1101
strumpack
Definition: StrumpackOptions.hpp:42
strumpack::SPOptions::use_METIS_NodeNDP
bool use_METIS_NodeNDP() const
Definition: StrumpackOptions.hpp:833
strumpack::MatchingJob::COMBBLAS
@ COMBBLAS
strumpack::MatchingJob::MAX_CARDINALITY
@ MAX_CARDINALITY
strumpack::SPOptions::set_GramSchmidt_type
void set_GramSchmidt_type(GramSchmidtType t)
Definition: StrumpackOptions.hpp:288
strumpack::KrylovSolver
KrylovSolver
Definition: StrumpackOptions.hpp:144
strumpack::default_abs_tol
real_t default_abs_tol()
Definition: StrumpackOptions.hpp:172
strumpack::SPOptions::enable_METIS_NodeNDP
void enable_METIS_NodeNDP()
Definition: StrumpackOptions.hpp:423
strumpack::SPOptions::HSS_options
HSS::HSSOptions< scalar_t > & HSS_options()
Definition: StrumpackOptions.hpp:1095
strumpack::CompressionType::LOSSY
@ LOSSY
strumpack::SPOptions::set_from_command_line
void set_from_command_line()
Definition: StrumpackOptions.hpp:1126
strumpack::CompressionType::HODLR
@ HODLR
strumpack::SPOptions
Options for the sparse solver.
Definition: StrumpackOptions.hpp:195
strumpack::SPOptions::set_reordering_method
void set_reordering_method(ReorderingStrategy m)
Definition: StrumpackOptions.hpp:302
strumpack::SPOptions::set_ny
void set_ny(int ny)
Definition: StrumpackOptions.hpp:367
strumpack::SPOptions::use_agg_amalg
bool use_agg_amalg() const
Definition: StrumpackOptions.hpp:852
strumpack::SPOptions::set_Krylov_solver
void set_Krylov_solver(KrylovSolver s)
Definition: StrumpackOptions.hpp:274
strumpack::CompressionType::BLR
@ BLR
strumpack::SPOptions::enable_agg_amalg
void enable_agg_amalg()
Definition: StrumpackOptions.hpp:482
strumpack::SPOptions< refine_t >::real_t
typename RealType< refine_t >::value_type real_t
Definition: StrumpackOptions.hpp:201
strumpack::SPOptions::compression_min_front_size
int compression_min_front_size(int l=0) const
Definition: StrumpackOptions.hpp:970
strumpack::SPOptions::set_dimensions
void set_dimensions(int nx, int ny=1, int nz=1)
Definition: StrumpackOptions.hpp:341
strumpack::SPOptions::set_maxit
void set_maxit(int maxit)
Definition: StrumpackOptions.hpp:238
strumpack::SPOptions::SPOptions
SPOptions()
Definition: StrumpackOptions.hpp:207
strumpack::SPOptions::set_rel_tol
void set_rel_tol(real_t rtol)
Definition: StrumpackOptions.hpp:246
strumpack::GramSchmidtType::MODIFIED
@ MODIFIED
strumpack::SPOptions::set_print_root_front_stats
void set_print_root_front_stats(bool b)
Definition: StrumpackOptions.hpp:734
strumpack::SPOptions::describe_options
void describe_options() const
strumpack::SPOptions::nd_param
int nd_param() const
Definition: StrumpackOptions.hpp:788
strumpack::SPOptions::nx
int nx() const
Definition: StrumpackOptions.hpp:802
strumpack::SPOptions::set_compression_rel_tol
void set_compression_rel_tol(real_t rtol)
Definition: StrumpackOptions.hpp:539
strumpack::HODLR::HODLROptions
Class containing several options for the HODLR code and data-structures.
Definition: HODLROptions.hpp:115
strumpack::CompressionType::BLR_HODLR
@ BLR_HODLR
strumpack::get_name
std::string get_name(ReorderingStrategy method)
strumpack::SPOptions::set_gmres_restart
void set_gmres_restart(int m)
Definition: StrumpackOptions.hpp:281
strumpack::get_description
std::string get_description(MatchingJob job)
strumpack::SPOptions::use_METIS_NodeND
bool use_METIS_NodeND() const
Definition: StrumpackOptions.hpp:839
strumpack::KrylovSolver::PREC_BICGSTAB
@ PREC_BICGSTAB
strumpack::SPOptions::compression_abs_tol
real_t compression_abs_tol(int l=0) const
Definition: StrumpackOptions.hpp:906
strumpack::SPOptions::maxit
int maxit() const
Definition: StrumpackOptions.hpp:746
strumpack::MatchingJob::NONE
@ NONE
strumpack::SPOptions::use_gpu
bool use_gpu() const
Definition: StrumpackOptions.hpp:1067
strumpack::ReorderingStrategy::RCM
@ RCM
strumpack::SPOptions::set_compression_leaf_size
void set_compression_leaf_size(int s)
Definition: StrumpackOptions.hpp:636
strumpack::SPOptions::set_matching
void set_matching(MatchingJob job)
Definition: StrumpackOptions.hpp:503
strumpack::SPOptions::enable_replace_tiny_pivots
void enable_replace_tiny_pivots()
Definition: StrumpackOptions.hpp:682
strumpack::SPOptions::set_compression_min_front_size
void set_compression_min_front_size(int s)
Definition: StrumpackOptions.hpp:607
strumpack::SPOptions::compression
CompressionType compression() const
Definition: StrumpackOptions.hpp:869
strumpack::MatchingJob::MAX_SMALLEST_DIAGONAL
@ MAX_SMALLEST_DIAGONAL
strumpack::SPOptions::gpu_streams
int gpu_streams() const
Definition: StrumpackOptions.hpp:1072
strumpack::SPOptions::separator_ordering_level
int separator_ordering_level() const
Definition: StrumpackOptions.hpp:1035
HSSOptions.hpp
Contains the HSSOptions class as well as general routines for HSS options.
strumpack::SPOptions::set_write_root_front
void set_write_root_front(bool b)
Definition: StrumpackOptions.hpp:707
strumpack::ReorderingStrategy::METIS
@ METIS
strumpack::SPOptions::compression_leaf_size
int compression_leaf_size(int l=0) const
Definition: StrumpackOptions.hpp:1010
strumpack::SPOptions::HSS_options
const HSS::HSSOptions< scalar_t > & HSS_options() const
Definition: StrumpackOptions.hpp:1089
strumpack::SPOptions::verbose
bool verbose() const
Definition: StrumpackOptions.hpp:740
strumpack::CompressionType::HSS
@ HSS
strumpack::SPOptions::disable_METIS_NodeND
void disable_METIS_NodeND()
Definition: StrumpackOptions.hpp:453
strumpack::SPOptions::replace_tiny_pivots
bool replace_tiny_pivots() const
Definition: StrumpackOptions.hpp:1047
strumpack::SPOptions::disable_gpu
void disable_gpu()
Definition: StrumpackOptions.hpp:718
strumpack::SPOptions::set_gpu_streams
void set_gpu_streams(int s)
Definition: StrumpackOptions.hpp:723
strumpack::SPOptions::enable_METIS_NodeND
void enable_METIS_NodeND()
Definition: StrumpackOptions.hpp:444
strumpack::ReorderingStrategy::PARMETIS
@ PARMETIS
strumpack::SPOptions::set_pivot_threshold
void set_pivot_threshold(real_t thresh)
Definition: StrumpackOptions.hpp:700
strumpack::ReorderingStrategy::SCOTCH
@ SCOTCH
strumpack::SPOptions::write_root_front
bool write_root_front() const
Definition: StrumpackOptions.hpp:1062
strumpack::SPOptions::separator_width
int separator_width() const
Definition: StrumpackOptions.hpp:827
strumpack::SPOptions::pivot_threshold
real_t pivot_threshold() const
Definition: StrumpackOptions.hpp:1057
strumpack::is_parallel
bool is_parallel(ReorderingStrategy method)
strumpack::SPOptions::set_compression_abs_tol
void set_compression_abs_tol(real_t atol)
Definition: StrumpackOptions.hpp:558
HODLROptions.hpp
Contains the class holding HODLR matrix options.
strumpack::SPOptions::gmres_restart
int gmres_restart() const
Definition: StrumpackOptions.hpp:770
strumpack::SPOptions::enable_gpu
void enable_gpu()
Definition: StrumpackOptions.hpp:713
strumpack::SPOptions::ny
int ny() const
Definition: StrumpackOptions.hpp:808
strumpack::ReorderingStrategy::GEOMETRIC
@ GEOMETRIC
strumpack::SPOptions::nd_planar_levels
int nd_planar_levels() const
Definition: StrumpackOptions.hpp:796
strumpack::SPOptions::SPOptions
SPOptions(int argc, const char *const argv[])
Definition: StrumpackOptions.hpp:220
strumpack::KrylovSolver::BICGSTAB
@ BICGSTAB
strumpack::SPOptions::print_root_front_stats
bool print_root_front_stats() const
Definition: StrumpackOptions.hpp:1083
strumpack::SPOptions::matching
MatchingJob matching() const
Definition: StrumpackOptions.hpp:858
strumpack::SPOptions::disable_assembly_tree_log
void disable_assembly_tree_log()
Definition: StrumpackOptions.hpp:514
strumpack::CompressionType
CompressionType
Definition: StrumpackOptions.hpp:74
strumpack::SPOptions::set_abs_tol
void set_abs_tol(real_t atol)
Definition: StrumpackOptions.hpp:257
strumpack::SPOptions::compression_min_sep_size
int compression_min_sep_size(int l=0) const
Definition: StrumpackOptions.hpp:931
strumpack::KrylovSolver::REFINE
@ REFINE
strumpack::ReorderingStrategy
ReorderingStrategy
Definition: StrumpackOptions.hpp:48
strumpack::SPOptions::set_compression
void set_compression(CompressionType c)
Definition: StrumpackOptions.hpp:524
strumpack::CompressionType::NONE
@ NONE
strumpack::MatchingJob
MatchingJob
Definition: StrumpackOptions.hpp:98