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 // #include "sparse/ordering/spectral/NDOptions.hpp"
42 
43 namespace strumpack {
44 
50  enum class ProportionalMapping {
51  FLOPS,
54  };
55 
60  enum class ReorderingStrategy {
61  NATURAL,
62  METIS,
63  PARMETIS,
64  SCOTCH,
65  PTSCOTCH,
66  RCM,
67  GEOMETRIC,
69  AMD,
70  MMD,
71  AND,
72  MLF,
73  SPECTRAL
74  };
75 
79  std::string get_name(ReorderingStrategy method);
80 
85 
91  enum class CompressionType {
92  NONE,
93  HSS,
94  BLR,
95  HODLR,
97  BLR_HODLR,
100  ZFP_BLR_HODLR,
104  LOSSLESS,
105  LOSSY
106  };
107 
111  std::string get_name(CompressionType comp);
112 
113 
119  enum class MatchingJob {
120  NONE,
128  COMBBLAS
129  };
130 
131  enum class EquilibrationType : char
132  { NONE='N', ROW='R', COLUMN='C', BOTH='B' };
133 
134 
139 
145 
149  std::string get_description(MatchingJob job);
150 
151 
156  enum class GramSchmidtType {
157  CLASSICAL,
158  MODIFIED
159  };
160 
165  enum class KrylovSolver {
166  AUTO,
168  DIRECT,
170  REFINE,
171  PREC_GMRES,
173  GMRES,
174  PREC_BICGSTAB,
176  BICGSTAB
177  };
178 
185  template<typename real_t> inline real_t default_rel_tol()
186  { return real_t(1.e-6); }
193  template<typename real_t> inline real_t default_abs_tol()
194  { return real_t(1.e-10); }
195  template<> inline float default_rel_tol() { return 1.e-4; }
196  template<> inline float default_abs_tol() { return 1.e-6; }
197 
198  inline int default_gpu_streams() { return 4; }
199 
216  template<typename scalar_t> class SPOptions {
217  public:
222  using real_t = typename RealType<scalar_t>::value_type;
223 
228  SPOptions() { hss_opts_.set_verbose(false); }
229 
241  SPOptions(int argc, const char* const argv[]) : argc_(argc), argv_(argv) {
242  hss_opts_.set_verbose(false);
243  blr_opts_.set_verbose(false);
244  hodlr_opts_.set_verbose(false);
245  // nd_opts_.set_verbose(false);
246  }
247 
254  void set_verbose(bool verbose) { verbose_ = verbose; }
255 
260  void set_maxit(int maxit) { assert(maxit >= 1); maxit_ = maxit; }
261 
268  void set_rel_tol(real_t rtol) {
269  assert(rtol <= real_t(1.) && rtol >= real_t(0.));
270  rel_tol_ = rtol;
271  }
272 
279  void set_abs_tol(real_t atol) {
280  assert(atol >= real_t(0.));
281  abs_tol_ = atol;
282  }
283 
296  void set_Krylov_solver(KrylovSolver s) { Krylov_solver_ = s; }
297 
303  void set_gmres_restart(int m) { assert(m >= 1); gmres_restart_ = m; }
304 
310  void set_GramSchmidt_type(GramSchmidtType t) { Gram_Schmidt_type_ = t; }
311 
324  void set_reordering_method(ReorderingStrategy m) { reordering_method_ = m; }
325 
333  { assert(nd_param>=0); nd_param_ = nd_param; }
334 
342  { assert(nd_planar_levels>=0); nd_planar_levels_ = nd_planar_levels; }
343 
363  void set_dimensions(int nx, int ny=1, int nz=1) {
364  assert(nx>=1 && ny>=1 && nz>=1);
365  nx_ = nx; ny_ = ny; nz_ = nz;
366  }
367 
377  void set_nx(int nx) {assert(nx>=1); nx_ = nx; }
378 
389  void set_ny(int ny) {assert(ny>=1); ny_ = ny; }
390 
401  void set_nz(int nz) { assert(nz>=1); nz_ = nz; }
402 
414  { assert(components>=1); components_ = components; }
415 
433  void set_separator_width(int width)
434  { assert(width>=1); separator_width_ = width; }
435 
445  void enable_METIS_NodeNDP() { use_METIS_NodeNDP_ = true; }
446 
456  void disable_METIS_NodeNDP() { use_METIS_NodeNDP_ = false; }
457 
458 
466  void enable_METIS_NodeND() { use_METIS_NodeNDP_ = false; }
467 
475  void disable_METIS_NodeND() { use_METIS_NodeNDP_ = true; }
476 
484  void enable_MUMPS_SYMQAMD() { use_MUMPS_SYMQAMD_ = true; }
485 
494  void disable_MUMPS_SYMQAMD() { use_MUMPS_SYMQAMD_ = false; }
495 
504  void enable_agg_amalg() { use_agg_amalg_ = true; }
505 
514  void disable_agg_amalg() { use_agg_amalg_ = false; }
515 
525  void set_matching(MatchingJob job) { matching_job_ = job; }
526 
530  void enable_assembly_tree_log() { log_assembly_tree_ = true; }
531 
536  void disable_assembly_tree_log() { log_assembly_tree_ = false; }
537 
546  void set_compression(CompressionType c) { comp_ = c; }
547 
562  hss_opts_.set_rel_tol(rtol);
563  blr_opts_.set_rel_tol(rtol);
564  hodlr_opts_.set_rel_tol(rtol);
565  }
566 
581  hss_opts_.set_abs_tol(atol);
582  blr_opts_.set_abs_tol(atol);
583  hodlr_opts_.set_abs_tol(atol);
584  }
585 
597  assert(s >= 0);
598  hss_min_sep_size_ = s;
599  blr_min_sep_size_ = s;
600  hodlr_min_sep_size_ = s;
601  lossy_min_sep_size_ = s;
602  }
603  void set_hss_min_sep_size(int s) {
604  assert(s >= 0);
605  hss_min_sep_size_ = s;
606  }
607  void set_hodlr_min_sep_size(int s) {
608  assert(s >= 0);
609  hodlr_min_sep_size_ = s;
610  }
611  void set_blr_min_sep_size(int s) {
612  assert(s >= 0);
613  blr_min_sep_size_ = s;
614  }
615  void set_lossy_min_sep_size(int s) {
616  assert(s >= 0);
617  lossy_min_sep_size_ = s;
618  }
619 
630  assert(s >= 0);
631  hss_min_front_size_ = s;
632  blr_min_front_size_ = s;
633  hodlr_min_front_size_ = s;
634  lossy_min_front_size_ = s;
635  }
636  void set_hss_min_front_size(int s) {
637  assert(s >= 0);
638  hss_min_front_size_ = s;
639  }
640  void set_hodlr_min_front_size(int s) {
641  assert(s >= 0);
642  hodlr_min_front_size_ = s;
643  }
644  void set_blr_min_front_size(int s) {
645  assert(s >= 0);
646  blr_min_front_size_ = s;
647  }
648  void set_lossy_min_front_size(int s) {
649  assert(s >= 0);
650  lossy_min_front_size_ = s;
651  }
652 
659  hss_opts_.set_leaf_size(s);
660  blr_opts_.set_leaf_size(s);
661  hodlr_opts_.set_leaf_size(s);
662  }
663 
680  { assert(l >= 0); sep_order_level_ = l; }
681 
682 #ifndef DOXYGEN_SHOULD_SKIP_THIS
683  void enable_indirect_sampling() { indirect_sampling_ = true; }
684  void disable_indirect_sampling() { indirect_sampling_ = false; }
685 #endif // DOXYGEN_SHOULD_SKIP_THIS
686 
704  void enable_replace_tiny_pivots() { replace_tiny_pivots_ = true; }
705 
711  void disable_replace_tiny_pivots() { replace_tiny_pivots_ = false; }
712 
722  void set_pivot_threshold(real_t thresh) { pivot_ = thresh; }
723 
729  void set_write_root_front(bool b) { write_root_front_ = b; }
730 
735  void enable_gpu() { use_gpu_ = true; }
736 
740  void disable_gpu() { use_gpu_ = false; }
741 
745  void set_gpu_streams(int s) { gpu_streams_ = s; }
746 
752  void enable_openmp_tree() { use_openmp_tree_ = true; }
753 
759  void disable_openmp_tree() { use_openmp_tree_ = false; }
760 
764  void set_lossy_precision(int p) { lossy_precision_ = p; }
765 
770  void set_print_compressed_front_stats(bool b) { print_comp_front_stats_ = b; }
771 
775  void set_proportional_mapping(ProportionalMapping pmap) { prop_map_ = pmap; }
776 
781  bool verbose() const { return verbose_; }
782 
787  int maxit() const { return maxit_; }
788 
793  real_t rel_tol() const { return rel_tol_; }
794 
799  real_t abs_tol() const { return abs_tol_; }
800 
805  KrylovSolver Krylov_solver() const { return Krylov_solver_; }
806 
811  int gmres_restart() const { return gmres_restart_; }
812 
817  GramSchmidtType GramSchmidt_type() const { return Gram_Schmidt_type_; }
818 
823  ReorderingStrategy reordering_method() const { return reordering_method_; }
824 
829  int nd_param() const { return nd_param_; }
830 
837  int nd_planar_levels() const { return nd_planar_levels_; }
838 
843  int nx() const { return nx_; }
844 
849  int ny() const { return ny_; }
850 
855  int nz() const { return nz_; }
856 
862  int components() const { return components_; }
863 
868  int separator_width() const { return separator_width_; }
869 
874  bool use_METIS_NodeNDP() const { return use_METIS_NodeNDP_; }
875 
880  bool use_METIS_NodeND() const { return !use_METIS_NodeNDP_; }
881 
886  bool use_MUMPS_SYMQAMD() const { return use_MUMPS_SYMQAMD_; }
887 
893  bool use_agg_amalg() const { return use_agg_amalg_; }
894 
899  MatchingJob matching() const { return matching_job_; }
900 
905  bool log_assembly_tree() const { return log_assembly_tree_; }
906 
910  CompressionType compression() const { return comp_; }
911 
912 
921  real_t compression_rel_tol(int l=0) const {
922  switch (comp_) {
924  return hss_opts_.rel_tol();
926  return blr_opts_.rel_tol();
928  return hodlr_opts_.rel_tol();
930  if (l==0) return hodlr_opts_.rel_tol();
931  else return blr_opts_.rel_tol();
933  if (l==0) return hodlr_opts_.rel_tol();
934  else return blr_opts_.rel_tol();
938  default: return 0.;
939  }
940  }
941 
950  real_t compression_abs_tol(int l=0) const {
951  switch (comp_) {
953  return hss_opts_.abs_tol();
955  return blr_opts_.abs_tol();
957  return hodlr_opts_.abs_tol();
959  if (l==0) return hodlr_opts_.abs_tol();
960  else return blr_opts_.abs_tol();
962  if (l==0) return hodlr_opts_.abs_tol();
963  else return blr_opts_.abs_tol();
967  default: return 0.;
968  }
969  }
970 
978  int compression_min_sep_size(int l=0) const {
979  switch (comp_) {
981  return hss_min_sep_size_;
983  return blr_min_sep_size_;
985  return hodlr_min_sep_size_;
987  if (l==0) return hodlr_min_sep_size_;
988  else return blr_min_sep_size_;
990  if (l==0) return hodlr_min_sep_size_;
991  else if (l==1) return blr_min_sep_size_;
992  else return lossy_min_sep_size_;
995  return lossy_min_sep_size_;
997  default:
998  return std::numeric_limits<int>::max();
999  }
1000  }
1001  int hss_min_sep_size() const {
1002  return hss_min_sep_size_;
1003  }
1004  int hodlr_min_sep_size() const {
1005  return hodlr_min_sep_size_;
1006  }
1007  int blr_min_sep_size() const {
1008  return blr_min_sep_size_;
1009  }
1010  int lossy_min_sep_size() const {
1011  return lossy_min_sep_size_;
1012  }
1013 
1021  int compression_min_front_size(int l=0) const {
1022  switch (comp_) {
1023  case CompressionType::HSS:
1024  return hss_min_front_size_;
1025  case CompressionType::BLR:
1026  return blr_min_front_size_;
1028  return hodlr_min_front_size_;
1030  if (l==0) return hodlr_min_front_size_;
1031  else return blr_min_front_size_;
1033  if (l==0) return hodlr_min_front_size_;
1034  else if (l==1) return blr_min_front_size_;
1035  else return lossy_min_front_size_;
1038  return lossy_min_front_size_;
1039  case CompressionType::NONE:
1040  default:
1041  return std::numeric_limits<int>::max();
1042  }
1043  }
1044  int hss_min_front_size() const {
1045  return hss_min_front_size_;
1046  }
1047  int hodlr_min_front_size() const {
1048  return hodlr_min_front_size_;
1049  }
1050  int blr_min_front_size() const {
1051  return blr_min_front_size_;
1052  }
1053  int lossy_min_front_size() const {
1054  return lossy_min_front_size_;
1055  }
1056 
1065  int compression_leaf_size(int l=0) const {
1066  switch (comp_) {
1067  case CompressionType::HSS:
1068  return hss_opts_.leaf_size();
1069  case CompressionType::BLR:
1070  return blr_opts_.leaf_size();
1072  return hodlr_opts_.leaf_size();
1074  if (l==0) return hodlr_opts_.leaf_size();
1075  else return blr_opts_.leaf_size();
1077  if (l==0) return hodlr_opts_.leaf_size();
1078  else if (l==1) return blr_opts_.leaf_size();
1079  else return 4;
1082  return 4;
1083  case CompressionType::NONE:
1084  default:
1085  return std::numeric_limits<int>::max();
1086  }
1087  }
1088 
1094  int separator_ordering_level() const { return sep_order_level_; }
1095 
1099  bool indirect_sampling() const { return indirect_sampling_; }
1100 
1106  bool replace_tiny_pivots() const { return replace_tiny_pivots_; }
1107 
1116  real_t pivot_threshold() const { return pivot_; }
1117 
1121  bool write_root_front() const { return write_root_front_; }
1122 
1126  bool use_gpu() const { return use_gpu_; }
1127 
1132  bool use_openmp_tree() const { return use_openmp_tree_; }
1133 
1137  int gpu_streams() const { return gpu_streams_; }
1138 
1142  int lossy_precision() const {
1143  return (compression() == CompressionType::LOSSLESS) ?
1144  -1 : lossy_precision_;
1145  }
1146 
1151  bool print_compressed_front_stats() const { return print_comp_front_stats_; }
1152 
1156  ProportionalMapping proportional_mapping() const { return prop_map_; }
1157 
1162  const HSS::HSSOptions<scalar_t>& HSS_options() const { return hss_opts_; }
1163 
1168  HSS::HSSOptions<scalar_t>& HSS_options() { return hss_opts_; }
1169 
1174  const BLR::BLROptions<scalar_t>& BLR_options() const { return blr_opts_; }
1175 
1180  BLR::BLROptions<scalar_t>& BLR_options() { return blr_opts_; }
1181 
1186  const HODLR::HODLROptions<scalar_t>& HODLR_options() const { return hodlr_opts_; }
1187 
1193 
1194  // /**
1195  // * Get a (const) reference to an object holding various options
1196  // * pertaining to the spectral nested dissection code.
1197  // */
1198  // const ordering::NDOptions& ND_options() const { return nd_opts_; }
1199 
1200  // /**
1201  // * Get a reference to an object holding various options pertaining
1202  // * to the spectral nested dissection code.
1203  // */
1204  // ordering::NDOptions& ND_options() { return nd_opts_; }
1205 
1212 
1225  void set_from_command_line(int argc, const char* const* cargv);
1226 
1231  void describe_options() const;
1232 
1233  private:
1234  bool verbose_ = true;
1236  int maxit_ = 5000;
1237  real_t rel_tol_ = default_rel_tol<real_t>();
1238  real_t abs_tol_ = default_abs_tol<real_t>();
1239  KrylovSolver Krylov_solver_ = KrylovSolver::AUTO;
1240  int gmres_restart_ = 30;
1241  GramSchmidtType Gram_Schmidt_type_ = GramSchmidtType::MODIFIED;
1243  ReorderingStrategy reordering_method_ = ReorderingStrategy::METIS;
1244  int nd_planar_levels_ = 0;
1245  int nd_param_ = 8;
1246  int nx_ = 1;
1247  int ny_ = 1;
1248  int nz_ = 1;
1249  int components_ = 1;
1250  int separator_width_ = 1;
1251  bool use_METIS_NodeNDP_ = false;
1252  bool use_MUMPS_SYMQAMD_ = false;
1253  bool use_agg_amalg_ = false;
1255  bool log_assembly_tree_ = false;
1256  bool replace_tiny_pivots_ = false;
1257  real_t pivot_ = std::sqrt(blas::lamch<real_t>('E'));
1258  bool write_root_front_ = false;
1259  bool print_comp_front_stats_ = false;
1261  bool use_openmp_tree_ = true;
1262 
1264 #if defined(STRUMPACK_USE_CUDA) || defined(STRUMPACK_USE_HIP) || defined(STRUMPACK_USE_SYCL)
1265  bool use_gpu_ = true;
1266 #else
1267  bool use_gpu_ = false;
1268 #endif
1269  int gpu_streams_ = default_gpu_streams();
1270 
1273 
1275  int hss_min_front_size_ = 100000;
1276  int hss_min_sep_size_ = 1000;
1277  int sep_order_level_ = 1;
1278  bool indirect_sampling_ = false;
1279  HSS::HSSOptions<scalar_t> hss_opts_;
1280 
1282  BLR::BLROptions<scalar_t> blr_opts_;
1283  int blr_min_front_size_ = 100000;
1284  int blr_min_sep_size_ = 512;
1285 
1287  HODLR::HODLROptions<scalar_t> hodlr_opts_;
1288  int hodlr_min_front_size_ = 100000;
1289  int hodlr_min_sep_size_ = 5000;
1290 
1292  int lossy_min_front_size_ = 100000;
1293  int lossy_min_sep_size_ = 8;
1294  int lossy_precision_ = 16;
1295 
1296  // ordering::NDOptions nd_opts_;
1297 
1298  int argc_ = 0;
1299  const char* const* argv_ = nullptr;
1300  };
1301 
1302 } // end namespace strumpack
1303 
1304 #endif // SPOPTIONS_HPP
Contains class holding BLROptions.
Contains the class holding HODLR matrix options.
Contains the HSSOptions class as well as general routines for HSS options.
Class containing several options for the BLR code and data-structures.
Definition: BLROptions.hpp:82
Class containing several options for the HODLR code and data-structures.
Definition: HODLROptions.hpp:117
Class containing several options for the HSS code and data-structures.
Definition: HSSOptions.hpp:152
Options for the sparse solver.
Definition: StrumpackOptions.hpp:216
void disable_assembly_tree_log()
Definition: StrumpackOptions.hpp:536
void enable_METIS_NodeND()
Definition: StrumpackOptions.hpp:466
int nd_planar_levels() const
Definition: StrumpackOptions.hpp:837
void set_from_command_line()
Definition: StrumpackOptions.hpp:1211
void set_gpu_streams(int s)
Definition: StrumpackOptions.hpp:745
int compression_min_front_size(int l=0) const
Definition: StrumpackOptions.hpp:1021
void set_matching(MatchingJob job)
Definition: StrumpackOptions.hpp:525
void set_verbose(bool verbose)
Definition: StrumpackOptions.hpp:254
int gmres_restart() const
Definition: StrumpackOptions.hpp:811
void disable_agg_amalg()
Definition: StrumpackOptions.hpp:514
int compression_min_sep_size(int l=0) const
Definition: StrumpackOptions.hpp:978
void set_reordering_method(ReorderingStrategy m)
Definition: StrumpackOptions.hpp:324
real_t pivot_threshold() const
Definition: StrumpackOptions.hpp:1116
void set_maxit(int maxit)
Definition: StrumpackOptions.hpp:260
void set_nd_planar_levels(int nd_planar_levels)
Definition: StrumpackOptions.hpp:341
int gpu_streams() const
Definition: StrumpackOptions.hpp:1137
void disable_openmp_tree()
Definition: StrumpackOptions.hpp:759
void set_compression_rel_tol(real_t rtol)
Definition: StrumpackOptions.hpp:561
void set_rel_tol(real_t rtol)
Definition: StrumpackOptions.hpp:268
void set_compression_leaf_size(int s)
Definition: StrumpackOptions.hpp:658
void disable_METIS_NodeNDP()
Definition: StrumpackOptions.hpp:456
void set_compression_min_sep_size(int s)
Definition: StrumpackOptions.hpp:596
bool verbose() const
Definition: StrumpackOptions.hpp:781
SPOptions()
Definition: StrumpackOptions.hpp:228
void disable_gpu()
Definition: StrumpackOptions.hpp:740
bool use_METIS_NodeND() const
Definition: StrumpackOptions.hpp:880
bool write_root_front() const
Definition: StrumpackOptions.hpp:1121
int separator_width() const
Definition: StrumpackOptions.hpp:868
typename RealType< scalar_t >::value_type real_t
Definition: StrumpackOptions.hpp:222
void set_nd_param(int nd_param)
Definition: StrumpackOptions.hpp:332
bool use_agg_amalg() const
Definition: StrumpackOptions.hpp:893
void set_lossy_precision(int p)
Definition: StrumpackOptions.hpp:764
void describe_options() const
CompressionType compression() const
Definition: StrumpackOptions.hpp:910
int nx() const
Definition: StrumpackOptions.hpp:843
void set_print_compressed_front_stats(bool b)
Definition: StrumpackOptions.hpp:770
void set_ny(int ny)
Definition: StrumpackOptions.hpp:389
void set_GramSchmidt_type(GramSchmidtType t)
Definition: StrumpackOptions.hpp:310
int lossy_precision() const
Definition: StrumpackOptions.hpp:1142
void set_from_command_line(int argc, const char *const *cargv)
void set_pivot_threshold(real_t thresh)
Definition: StrumpackOptions.hpp:722
void set_separator_width(int width)
Definition: StrumpackOptions.hpp:433
void set_abs_tol(real_t atol)
Definition: StrumpackOptions.hpp:279
void set_write_root_front(bool b)
Definition: StrumpackOptions.hpp:729
void disable_replace_tiny_pivots()
Definition: StrumpackOptions.hpp:711
void enable_openmp_tree()
Definition: StrumpackOptions.hpp:752
int ny() const
Definition: StrumpackOptions.hpp:849
bool use_gpu() const
Definition: StrumpackOptions.hpp:1126
void enable_gpu()
Definition: StrumpackOptions.hpp:735
void set_separator_ordering_level(int l)
Definition: StrumpackOptions.hpp:679
int nd_param() const
Definition: StrumpackOptions.hpp:829
bool replace_tiny_pivots() const
Definition: StrumpackOptions.hpp:1106
MatchingJob matching() const
Definition: StrumpackOptions.hpp:899
void enable_replace_tiny_pivots()
Definition: StrumpackOptions.hpp:704
void set_nz(int nz)
Definition: StrumpackOptions.hpp:401
real_t compression_rel_tol(int l=0) const
Definition: StrumpackOptions.hpp:921
ReorderingStrategy reordering_method() const
Definition: StrumpackOptions.hpp:823
bool indirect_sampling() const
Definition: StrumpackOptions.hpp:1099
int maxit() const
Definition: StrumpackOptions.hpp:787
bool log_assembly_tree() const
Definition: StrumpackOptions.hpp:905
void enable_agg_amalg()
Definition: StrumpackOptions.hpp:504
int compression_leaf_size(int l=0) const
Definition: StrumpackOptions.hpp:1065
void disable_METIS_NodeND()
Definition: StrumpackOptions.hpp:475
void set_gmres_restart(int m)
Definition: StrumpackOptions.hpp:303
bool use_METIS_NodeNDP() const
Definition: StrumpackOptions.hpp:874
bool use_openmp_tree() const
Definition: StrumpackOptions.hpp:1132
const HODLR::HODLROptions< scalar_t > & HODLR_options() const
Definition: StrumpackOptions.hpp:1186
void disable_MUMPS_SYMQAMD()
Definition: StrumpackOptions.hpp:494
int components() const
Definition: StrumpackOptions.hpp:862
BLR::BLROptions< scalar_t > & BLR_options()
Definition: StrumpackOptions.hpp:1180
void set_compression_abs_tol(real_t atol)
Definition: StrumpackOptions.hpp:580
KrylovSolver Krylov_solver() const
Definition: StrumpackOptions.hpp:805
void set_dimensions(int nx, int ny=1, int nz=1)
Definition: StrumpackOptions.hpp:363
HODLR::HODLROptions< scalar_t > & HODLR_options()
Definition: StrumpackOptions.hpp:1192
int nz() const
Definition: StrumpackOptions.hpp:855
void set_nx(int nx)
Definition: StrumpackOptions.hpp:377
const BLR::BLROptions< scalar_t > & BLR_options() const
Definition: StrumpackOptions.hpp:1174
const HSS::HSSOptions< scalar_t > & HSS_options() const
Definition: StrumpackOptions.hpp:1162
HSS::HSSOptions< scalar_t > & HSS_options()
Definition: StrumpackOptions.hpp:1168
void set_compression(CompressionType c)
Definition: StrumpackOptions.hpp:546
bool print_compressed_front_stats() const
Definition: StrumpackOptions.hpp:1151
int separator_ordering_level() const
Definition: StrumpackOptions.hpp:1094
real_t abs_tol() const
Definition: StrumpackOptions.hpp:799
real_t rel_tol() const
Definition: StrumpackOptions.hpp:793
void set_Krylov_solver(KrylovSolver s)
Definition: StrumpackOptions.hpp:296
void enable_METIS_NodeNDP()
Definition: StrumpackOptions.hpp:445
GramSchmidtType GramSchmidt_type() const
Definition: StrumpackOptions.hpp:817
void set_compression_min_front_size(int s)
Definition: StrumpackOptions.hpp:629
void enable_assembly_tree_log()
Definition: StrumpackOptions.hpp:530
void enable_MUMPS_SYMQAMD()
Definition: StrumpackOptions.hpp:484
real_t compression_abs_tol(int l=0) const
Definition: StrumpackOptions.hpp:950
SPOptions(int argc, const char *const argv[])
Definition: StrumpackOptions.hpp:241
bool use_MUMPS_SYMQAMD() const
Definition: StrumpackOptions.hpp:886
ProportionalMapping proportional_mapping() const
Definition: StrumpackOptions.hpp:1156
void set_proportional_mapping(ProportionalMapping pmap)
Definition: StrumpackOptions.hpp:775
void set_components(int components)
Definition: StrumpackOptions.hpp:413
Definition: StrumpackOptions.hpp:43
GramSchmidtType
Definition: StrumpackOptions.hpp:156
real_t default_rel_tol()
Definition: StrumpackOptions.hpp:185
std::string get_description(MatchingJob job)
bool is_parallel(ReorderingStrategy method)
ProportionalMapping
Definition: StrumpackOptions.hpp:50
MatchingJob
Definition: StrumpackOptions.hpp:119
CompressionType
Definition: StrumpackOptions.hpp:91
ReorderingStrategy
Definition: StrumpackOptions.hpp:60
KrylovSolver
Definition: StrumpackOptions.hpp:165
std::string get_name(ReorderingStrategy method)
MatchingJob get_matching(int job)
real_t default_abs_tol()
Definition: StrumpackOptions.hpp:193