Loading...
Searching...
No Matches
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#include <cstdlib>
37
38#include "dense/BLASLAPACKWrapper.hpp"
39#include "HSS/HSSOptions.hpp"
40#include "BLR/BLROptions.hpp"
42// #include "sparse/ordering/spectral/NDOptions.hpp"
43
44namespace strumpack {
45
52 FLOPS,
55 };
56
61 enum class ReorderingStrategy {
62 NATURAL,
63 METIS,
64 PARMETIS,
65 SCOTCH,
66 PTSCOTCH,
67 RCM,
68 GEOMETRIC,
70 AMD,
71 MMD,
72 AND,
73 MLF,
75 };
76
80 std::string get_name(ReorderingStrategy method);
81
86
92 enum class CompressionType {
93 NONE,
94 HSS,
95 BLR,
96 HODLR,
98 BLR_HODLR,
105 LOSSLESS,
106 LOSSY
107 };
108
112 std::string get_name(CompressionType comp);
113
114
131
132 enum class EquilibrationType : char
133 { NONE='N', ROW='R', COLUMN='C', BOTH='B' };
134
135
140
146
150 std::string get_description(MatchingJob job);
151
152
157 enum class GramSchmidtType {
158 CLASSICAL,
159 MODIFIED
160 };
161
166 enum class KrylovSolver {
167 AUTO,
169 DIRECT,
171 REFINE,
172 PREC_GMRES,
174 GMRES,
177 BICGSTAB
178 };
179
186 template<typename real_t> inline real_t default_rel_tol()
187 { return real_t(1.e-6); }
194 template<typename real_t> inline real_t default_abs_tol()
195 { return real_t(1.e-10); }
196 template<> inline float default_rel_tol() { return 1.e-4; }
197 template<> inline float default_abs_tol() { return 1.e-6; }
198
199 inline int default_gpu_streams() { return 4; }
200
217 template<typename scalar_t> class SPOptions {
218 public:
223 using real_t = typename RealType<scalar_t>::value_type;
224
230 hss_opts_.set_verbose(false);
231 blr_opts_.set_verbose(false);
232 hodlr_opts_.set_verbose(false);
233 }
234
246 SPOptions(int argc, const char* const argv[]) : argc_(argc), argv_(argv) {
247 hss_opts_.set_verbose(false);
248 blr_opts_.set_verbose(false);
249 hodlr_opts_.set_verbose(false);
250 // nd_opts_.set_verbose(false);
251 }
252
259 void set_verbose(bool verbose) { verbose_ = verbose; }
260
265 void set_maxit(int maxit) { assert(maxit >= 1); maxit_ = maxit; }
266
273 void set_rel_tol(real_t rtol) {
274 assert(rtol <= real_t(1.) && rtol >= real_t(0.));
275 rel_tol_ = rtol;
276 }
277
284 void set_abs_tol(real_t atol) {
285 assert(atol >= real_t(0.));
286 abs_tol_ = atol;
287 }
288
301 void set_Krylov_solver(KrylovSolver s) { Krylov_solver_ = s; }
302
308 void set_gmres_restart(int m) { assert(m >= 1); gmres_restart_ = m; }
309
315 void set_GramSchmidt_type(GramSchmidtType t) { Gram_Schmidt_type_ = t; }
316
329 void set_reordering_method(ReorderingStrategy m) { reordering_method_ = m; }
330
338 { assert(nd_param>=0); nd_param_ = nd_param; }
339
347 { assert(nd_planar_levels>=0); nd_planar_levels_ = nd_planar_levels; }
348
368 void set_dimensions(int nx, int ny=1, int nz=1) {
369 assert(nx>=1 && ny>=1 && nz>=1);
370 nx_ = nx; ny_ = ny; nz_ = nz;
371 }
372
382 void set_nx(int nx) {assert(nx>=1); nx_ = nx; }
383
394 void set_ny(int ny) {assert(ny>=1); ny_ = ny; }
395
406 void set_nz(int nz) { assert(nz>=1); nz_ = nz; }
407
419 { assert(components>=1); components_ = components; }
420
438 void set_separator_width(int width)
439 { assert(width>=1); separator_width_ = width; }
440
450 void enable_METIS_NodeNDP() { use_METIS_NodeNDP_ = true; }
451
461 void disable_METIS_NodeNDP() { use_METIS_NodeNDP_ = false; }
462
463
471 void enable_METIS_NodeND() { use_METIS_NodeNDP_ = false; }
472
480 void disable_METIS_NodeND() { use_METIS_NodeNDP_ = true; }
481
489 void enable_MUMPS_SYMQAMD() { use_MUMPS_SYMQAMD_ = true; }
490
499 void disable_MUMPS_SYMQAMD() { use_MUMPS_SYMQAMD_ = false; }
500
509 void enable_agg_amalg() { use_agg_amalg_ = true; }
510
519 void disable_agg_amalg() { use_agg_amalg_ = false; }
520
530 void set_matching(MatchingJob job) { matching_job_ = job; }
531
535 void enable_assembly_tree_log() { log_assembly_tree_ = true; }
536
541 void disable_assembly_tree_log() { log_assembly_tree_ = false; }
542
551 void set_compression(CompressionType c) { comp_ = c; }
552
567 hss_opts_.set_rel_tol(rtol);
568 blr_opts_.set_rel_tol(rtol);
569 hodlr_opts_.set_rel_tol(rtol);
570 }
571
586 hss_opts_.set_abs_tol(atol);
587 blr_opts_.set_abs_tol(atol);
588 hodlr_opts_.set_abs_tol(atol);
589 }
590
602 assert(s >= 0);
603 hss_min_sep_size_ = s;
604 blr_min_sep_size_ = s;
605 hodlr_min_sep_size_ = s;
606 lossy_min_sep_size_ = s;
607 }
608 void set_hss_min_sep_size(int s) {
609 assert(s >= 0);
610 hss_min_sep_size_ = s;
611 }
612 void set_hodlr_min_sep_size(int s) {
613 assert(s >= 0);
614 hodlr_min_sep_size_ = s;
615 }
616 void set_blr_min_sep_size(int s) {
617 assert(s >= 0);
618 blr_min_sep_size_ = s;
619 }
620 void set_lossy_min_sep_size(int s) {
621 assert(s >= 0);
622 lossy_min_sep_size_ = s;
623 }
624
635 assert(s >= 0);
636 hss_min_front_size_ = s;
637 blr_min_front_size_ = s;
638 hodlr_min_front_size_ = s;
639 lossy_min_front_size_ = s;
640 }
641 void set_hss_min_front_size(int s) {
642 assert(s >= 0);
643 hss_min_front_size_ = s;
644 }
645 void set_hodlr_min_front_size(int s) {
646 assert(s >= 0);
647 hodlr_min_front_size_ = s;
648 }
649 void set_blr_min_front_size(int s) {
650 assert(s >= 0);
651 blr_min_front_size_ = s;
652 }
653 void set_lossy_min_front_size(int s) {
654 assert(s >= 0);
655 lossy_min_front_size_ = s;
656 }
657
664 hss_opts_.set_leaf_size(s);
665 blr_opts_.set_leaf_size(s);
666 hodlr_opts_.set_leaf_size(s);
667 }
668
685 { assert(l >= 0); sep_order_level_ = l; }
686
687#ifndef DOXYGEN_SHOULD_SKIP_THIS
688 void enable_indirect_sampling() { indirect_sampling_ = true; }
689 void disable_indirect_sampling() { indirect_sampling_ = false; }
690#endif // DOXYGEN_SHOULD_SKIP_THIS
691
709 void enable_replace_tiny_pivots() { replace_tiny_pivots_ = true; }
710
716 void disable_replace_tiny_pivots() { replace_tiny_pivots_ = false; }
717
728 pivot_ = thresh;
729 hss_opts_.set_pivot_threshold(thresh);
730 blr_opts_.set_pivot_threshold(thresh);
731 hodlr_opts_.set_pivot_threshold(thresh);
732 }
733
739 void set_write_root_front(bool b) { write_root_front_ = b; }
740
745 void enable_gpu() { use_gpu_ = true; }
746
750 void disable_gpu() { use_gpu_ = false; }
751
752// /**
753// * Disable symmetric solver.
754// */
755// void disable_symmetric() { use_symmetric_ = false; }
756//
757// /**
758// * Disable positive_definite solver.
759// */
760// void disable_positive_definite() { use_positive_definite_ = false; }
764 void enable_symmetric() { use_symmetric_ = true; }
765
766
770 void enable_positive_definite() { use_positive_definite_ = true; }
771
772
776 void set_gpu_streams(int s) { gpu_streams_ = s; }
777
783 void enable_openmp_tree() { use_openmp_tree_ = true; }
784
790 void disable_openmp_tree() { use_openmp_tree_ = false; }
791
801 lossy_precision_ = p;
802 }
803
809 void set_lossy_accuracy(double a) { lossy_accuracy_ = a; }
810
815 void set_print_compressed_front_stats(bool b) { print_comp_front_stats_ = b; }
816
820 void set_proportional_mapping(ProportionalMapping pmap) { prop_map_ = pmap; }
821
826 bool verbose() const { return verbose_; }
827
832 int maxit() const { return maxit_; }
833
838 real_t rel_tol() const { return rel_tol_; }
839
844 real_t abs_tol() const { return abs_tol_; }
845
850 KrylovSolver Krylov_solver() const { return Krylov_solver_; }
851
856 int gmres_restart() const { return gmres_restart_; }
857
862 GramSchmidtType GramSchmidt_type() const { return Gram_Schmidt_type_; }
863
868 ReorderingStrategy reordering_method() const { return reordering_method_; }
869
874 int nd_param() const { return nd_param_; }
875
882 int nd_planar_levels() const { return nd_planar_levels_; }
883
888 int nx() const { return nx_; }
889
894 int ny() const { return ny_; }
895
900 int nz() const { return nz_; }
901
907 int components() const { return components_; }
908
913 int separator_width() const { return separator_width_; }
914
919 bool use_METIS_NodeNDP() const { return use_METIS_NodeNDP_; }
920
925 bool use_METIS_NodeND() const { return !use_METIS_NodeNDP_; }
926
931 bool use_MUMPS_SYMQAMD() const { return use_MUMPS_SYMQAMD_; }
932
938 bool use_agg_amalg() const { return use_agg_amalg_; }
939
944 MatchingJob matching() const { return matching_job_; }
945
950 bool log_assembly_tree() const { return log_assembly_tree_; }
951
955 CompressionType compression() const { return comp_; }
956
957
967 switch (comp_) {
969 return hss_opts_.rel_tol();
971 return blr_opts_.rel_tol();
973 return hodlr_opts_.rel_tol();
975 if (l==0) return hodlr_opts_.rel_tol();
976 else return blr_opts_.rel_tol();
978 if (l==0) return hodlr_opts_.rel_tol();
979 else return blr_opts_.rel_tol();
983 default: return 0.;
984 }
985 }
986
996 switch (comp_) {
998 return hss_opts_.abs_tol();
1000 return blr_opts_.abs_tol();
1002 return hodlr_opts_.abs_tol();
1004 if (l==0) return hodlr_opts_.abs_tol();
1005 else return blr_opts_.abs_tol();
1007 if (l==0) return hodlr_opts_.abs_tol();
1008 else return blr_opts_.abs_tol();
1012 default: return 0.;
1013 }
1014 }
1015
1023 int compression_min_sep_size(int l=0) const {
1024 switch (comp_) {
1026 return hss_min_sep_size_;
1028 return blr_min_sep_size_;
1030 return hodlr_min_sep_size_;
1032 if (l==0) return hodlr_min_sep_size_;
1033 else return blr_min_sep_size_;
1035 if (l==0) return hodlr_min_sep_size_;
1036 else if (l==1) return blr_min_sep_size_;
1037 else return lossy_min_sep_size_;
1040 return lossy_min_sep_size_;
1042 default:
1043 return std::numeric_limits<int>::max();
1044 }
1045 }
1046 int hss_min_sep_size() const {
1047 return hss_min_sep_size_;
1048 }
1049 int hodlr_min_sep_size() const {
1050 return hodlr_min_sep_size_;
1051 }
1052 int blr_min_sep_size() const {
1053 return blr_min_sep_size_;
1054 }
1055 int lossy_min_sep_size() const {
1056 return lossy_min_sep_size_;
1057 }
1058
1066 int compression_min_front_size(int l=0) const {
1067 switch (comp_) {
1069 return hss_min_front_size_;
1071 return blr_min_front_size_;
1073 return hodlr_min_front_size_;
1075 if (l==0) return hodlr_min_front_size_;
1076 else return blr_min_front_size_;
1078 if (l==0) return hodlr_min_front_size_;
1079 else if (l==1) return blr_min_front_size_;
1080 else return lossy_min_front_size_;
1083 return lossy_min_front_size_;
1085 default:
1086 return std::numeric_limits<int>::max();
1087 }
1088 }
1089 int hss_min_front_size() const {
1090 return hss_min_front_size_;
1091 }
1092 int hodlr_min_front_size() const {
1093 return hodlr_min_front_size_;
1094 }
1095 int blr_min_front_size() const {
1096 return blr_min_front_size_;
1097 }
1098 int lossy_min_front_size() const {
1099 return lossy_min_front_size_;
1100 }
1101
1110 int compression_leaf_size(int l=0) const {
1111 switch (comp_) {
1113 return hss_opts_.leaf_size();
1115 return blr_opts_.leaf_size();
1117 return hodlr_opts_.leaf_size();
1119 if (l==0) return hodlr_opts_.leaf_size();
1120 else return blr_opts_.leaf_size();
1122 if (l==0) return hodlr_opts_.leaf_size();
1123 else if (l==1) return blr_opts_.leaf_size();
1124 else return 4;
1127 return 4;
1129 default:
1130 return std::numeric_limits<int>::max();
1131 }
1132 }
1133
1139 int separator_ordering_level() const { return sep_order_level_; }
1140
1144 bool indirect_sampling() const { return indirect_sampling_; }
1145
1151 bool replace_tiny_pivots() const { return replace_tiny_pivots_; }
1152
1161 real_t pivot_threshold() const { return pivot_; }
1162
1166 bool write_root_front() const { return write_root_front_; }
1167
1171 bool use_gpu() const { return use_gpu_; }
1172
1176 bool use_symmetric() const { return use_symmetric_; }
1177
1181 bool use_positive_definite() const { return use_positive_definite_; }
1182
1187 bool use_openmp_tree() const { return use_openmp_tree_; }
1188
1192 int gpu_streams() const { return gpu_streams_; }
1193
1197 bool use_gpu_aware_mpi() const { return use_gpu_aware_mpi_; }
1198
1202 int lossy_precision() const {
1204 -1 : lossy_precision_;
1205 }
1206
1210 double lossy_accuracy() const {
1212 -1 : lossy_accuracy_;
1213 }
1214
1219 bool print_compressed_front_stats() const { return print_comp_front_stats_; }
1220
1224 ProportionalMapping proportional_mapping() const { return prop_map_; }
1225
1230 const HSS::HSSOptions<scalar_t>& HSS_options() const { return hss_opts_; }
1231
1237
1242 const BLR::BLROptions<scalar_t>& BLR_options() const { return blr_opts_; }
1243
1249
1254 const HODLR::HODLROptions<scalar_t>& HODLR_options() const { return hodlr_opts_; }
1255
1261
1262 // /**
1263 // * Get a (const) reference to an object holding various options
1264 // * pertaining to the spectral nested dissection code.
1265 // */
1266 // const ordering::NDOptions& ND_options() const { return nd_opts_; }
1267
1268 // /**
1269 // * Get a reference to an object holding various options pertaining
1270 // * to the spectral nested dissection code.
1271 // */
1272 // ordering::NDOptions& ND_options() { return nd_opts_; }
1273
1280
1293 void set_from_command_line(int argc, const char* const* cargv);
1294
1299 void describe_options() const;
1300
1301 private:
1302 bool verbose_ = true;
1304 int maxit_ = 5000;
1305 real_t rel_tol_ = default_rel_tol<real_t>();
1306 real_t abs_tol_ = default_abs_tol<real_t>();
1307 KrylovSolver Krylov_solver_ = KrylovSolver::AUTO;
1308 int gmres_restart_ = 30;
1309 GramSchmidtType Gram_Schmidt_type_ = GramSchmidtType::MODIFIED;
1311 ReorderingStrategy reordering_method_ = ReorderingStrategy::METIS;
1312 int nd_planar_levels_ = 0;
1313 int nd_param_ = 8;
1314 int nx_ = 1;
1315 int ny_ = 1;
1316 int nz_ = 1;
1317 int components_ = 1;
1318 int separator_width_ = 1;
1319 bool use_METIS_NodeNDP_ = false;
1320 bool use_MUMPS_SYMQAMD_ = false;
1321 bool use_agg_amalg_ = false;
1323 bool log_assembly_tree_ = false;
1324 bool replace_tiny_pivots_ = false;
1325 real_t pivot_ = std::sqrt(blas::lamch<real_t>('E'));
1326 bool write_root_front_ = false;
1327 bool print_comp_front_stats_ = false;
1329 bool use_openmp_tree_ = true;
1330 bool use_symmetric_ = false;
1331 bool use_positive_definite_ = false;
1332
1334#if defined(STRUMPACK_USE_GPU)
1335 bool use_gpu_ = true;
1336#else
1337 bool use_gpu_ = false;
1338#endif
1339 bool use_gpu_aware_mpi_ = std::getenv("STRUMPACK_GPU_AWARE_MPI");
1340 int gpu_streams_ = default_gpu_streams();
1341
1344
1346 int hss_min_front_size_ = 100000;
1347 int hss_min_sep_size_ = 1000;
1348 int sep_order_level_ = 1;
1349 bool indirect_sampling_ = false;
1350 HSS::HSSOptions<scalar_t> hss_opts_;
1351
1353 BLR::BLROptions<scalar_t> blr_opts_;
1354 int blr_min_front_size_ = 100000;
1355 int blr_min_sep_size_ = 512;
1356
1359 int hodlr_min_front_size_ = 100000;
1360 int hodlr_min_sep_size_ = 5000;
1361
1363 int lossy_min_front_size_ = 100000;
1364 int lossy_min_sep_size_ = 8;
1365 int lossy_precision_ = 16;
1366 double lossy_accuracy_ = 1e-3;
1367
1368 // ordering::NDOptions nd_opts_;
1369
1370 int argc_ = 0;
1371 const char* const* argv_ = nullptr;
1372 };
1373
1374} // end namespace strumpack
1375
1376#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.
Definition BLRMatrix.hpp:50
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:217
void disable_assembly_tree_log()
Definition StrumpackOptions.hpp:541
void enable_METIS_NodeND()
Definition StrumpackOptions.hpp:471
double lossy_accuracy() const
Definition StrumpackOptions.hpp:1210
int nd_planar_levels() const
Definition StrumpackOptions.hpp:882
void set_from_command_line()
Definition StrumpackOptions.hpp:1279
void set_gpu_streams(int s)
Definition StrumpackOptions.hpp:776
int compression_min_front_size(int l=0) const
Definition StrumpackOptions.hpp:1066
void set_matching(MatchingJob job)
Definition StrumpackOptions.hpp:530
void set_verbose(bool verbose)
Definition StrumpackOptions.hpp:259
int gmres_restart() const
Definition StrumpackOptions.hpp:856
void disable_agg_amalg()
Definition StrumpackOptions.hpp:519
int compression_min_sep_size(int l=0) const
Definition StrumpackOptions.hpp:1023
void set_reordering_method(ReorderingStrategy m)
Definition StrumpackOptions.hpp:329
real_t pivot_threshold() const
Definition StrumpackOptions.hpp:1161
void set_maxit(int maxit)
Definition StrumpackOptions.hpp:265
void enable_symmetric()
Definition StrumpackOptions.hpp:764
HSS::HSSOptions< scalar_t > & HSS_options()
Definition StrumpackOptions.hpp:1236
void set_nd_planar_levels(int nd_planar_levels)
Definition StrumpackOptions.hpp:346
int gpu_streams() const
Definition StrumpackOptions.hpp:1192
void disable_openmp_tree()
Definition StrumpackOptions.hpp:790
void set_compression_rel_tol(real_t rtol)
Definition StrumpackOptions.hpp:566
void set_rel_tol(real_t rtol)
Definition StrumpackOptions.hpp:273
void set_compression_leaf_size(int s)
Definition StrumpackOptions.hpp:663
void disable_METIS_NodeNDP()
Definition StrumpackOptions.hpp:461
void set_compression_min_sep_size(int s)
Definition StrumpackOptions.hpp:601
bool verbose() const
Definition StrumpackOptions.hpp:826
SPOptions()
Definition StrumpackOptions.hpp:229
void disable_gpu()
Definition StrumpackOptions.hpp:750
void enable_positive_definite()
Definition StrumpackOptions.hpp:770
bool use_METIS_NodeND() const
Definition StrumpackOptions.hpp:925
bool write_root_front() const
Definition StrumpackOptions.hpp:1166
int separator_width() const
Definition StrumpackOptions.hpp:913
typename RealType< scalar_t >::value_type real_t
Definition StrumpackOptions.hpp:223
void set_nd_param(int nd_param)
Definition StrumpackOptions.hpp:337
bool use_agg_amalg() const
Definition StrumpackOptions.hpp:938
void set_lossy_precision(int p)
Definition StrumpackOptions.hpp:800
void describe_options() const
CompressionType compression() const
Definition StrumpackOptions.hpp:955
int nx() const
Definition StrumpackOptions.hpp:888
void set_print_compressed_front_stats(bool b)
Definition StrumpackOptions.hpp:815
void set_ny(int ny)
Definition StrumpackOptions.hpp:394
BLR::BLROptions< scalar_t > & BLR_options()
Definition StrumpackOptions.hpp:1248
void set_GramSchmidt_type(GramSchmidtType t)
Definition StrumpackOptions.hpp:315
int lossy_precision() const
Definition StrumpackOptions.hpp:1202
void set_from_command_line(int argc, const char *const *cargv)
void set_pivot_threshold(real_t thresh)
Definition StrumpackOptions.hpp:727
void set_separator_width(int width)
Definition StrumpackOptions.hpp:438
void set_abs_tol(real_t atol)
Definition StrumpackOptions.hpp:284
void set_write_root_front(bool b)
Definition StrumpackOptions.hpp:739
void disable_replace_tiny_pivots()
Definition StrumpackOptions.hpp:716
void enable_openmp_tree()
Definition StrumpackOptions.hpp:783
int ny() const
Definition StrumpackOptions.hpp:894
bool use_gpu() const
Definition StrumpackOptions.hpp:1171
void enable_gpu()
Definition StrumpackOptions.hpp:745
void set_separator_ordering_level(int l)
Definition StrumpackOptions.hpp:684
int nd_param() const
Definition StrumpackOptions.hpp:874
bool replace_tiny_pivots() const
Definition StrumpackOptions.hpp:1151
const BLR::BLROptions< scalar_t > & BLR_options() const
Definition StrumpackOptions.hpp:1242
MatchingJob matching() const
Definition StrumpackOptions.hpp:944
void enable_replace_tiny_pivots()
Definition StrumpackOptions.hpp:709
void set_nz(int nz)
Definition StrumpackOptions.hpp:406
real_t compression_rel_tol(int l=0) const
Definition StrumpackOptions.hpp:966
ReorderingStrategy reordering_method() const
Definition StrumpackOptions.hpp:868
bool indirect_sampling() const
Definition StrumpackOptions.hpp:1144
int maxit() const
Definition StrumpackOptions.hpp:832
void set_lossy_accuracy(double a)
Definition StrumpackOptions.hpp:809
bool log_assembly_tree() const
Definition StrumpackOptions.hpp:950
void enable_agg_amalg()
Definition StrumpackOptions.hpp:509
int compression_leaf_size(int l=0) const
Definition StrumpackOptions.hpp:1110
void disable_METIS_NodeND()
Definition StrumpackOptions.hpp:480
void set_gmres_restart(int m)
Definition StrumpackOptions.hpp:308
const HODLR::HODLROptions< scalar_t > & HODLR_options() const
Definition StrumpackOptions.hpp:1254
bool use_positive_definite() const
Definition StrumpackOptions.hpp:1181
bool use_METIS_NodeNDP() const
Definition StrumpackOptions.hpp:919
bool use_openmp_tree() const
Definition StrumpackOptions.hpp:1187
void disable_MUMPS_SYMQAMD()
Definition StrumpackOptions.hpp:499
int components() const
Definition StrumpackOptions.hpp:907
void set_compression_abs_tol(real_t atol)
Definition StrumpackOptions.hpp:585
KrylovSolver Krylov_solver() const
Definition StrumpackOptions.hpp:850
void set_dimensions(int nx, int ny=1, int nz=1)
Definition StrumpackOptions.hpp:368
bool use_symmetric() const
Definition StrumpackOptions.hpp:1176
int nz() const
Definition StrumpackOptions.hpp:900
void set_nx(int nx)
Definition StrumpackOptions.hpp:382
void set_compression(CompressionType c)
Definition StrumpackOptions.hpp:551
bool print_compressed_front_stats() const
Definition StrumpackOptions.hpp:1219
int separator_ordering_level() const
Definition StrumpackOptions.hpp:1139
real_t abs_tol() const
Definition StrumpackOptions.hpp:844
real_t rel_tol() const
Definition StrumpackOptions.hpp:838
bool use_gpu_aware_mpi() const
Definition StrumpackOptions.hpp:1197
void set_Krylov_solver(KrylovSolver s)
Definition StrumpackOptions.hpp:301
void enable_METIS_NodeNDP()
Definition StrumpackOptions.hpp:450
const HSS::HSSOptions< scalar_t > & HSS_options() const
Definition StrumpackOptions.hpp:1230
GramSchmidtType GramSchmidt_type() const
Definition StrumpackOptions.hpp:862
void set_compression_min_front_size(int s)
Definition StrumpackOptions.hpp:634
void enable_assembly_tree_log()
Definition StrumpackOptions.hpp:535
void enable_MUMPS_SYMQAMD()
Definition StrumpackOptions.hpp:489
real_t compression_abs_tol(int l=0) const
Definition StrumpackOptions.hpp:995
SPOptions(int argc, const char *const argv[])
Definition StrumpackOptions.hpp:246
bool use_MUMPS_SYMQAMD() const
Definition StrumpackOptions.hpp:931
ProportionalMapping proportional_mapping() const
Definition StrumpackOptions.hpp:1224
HODLR::HODLROptions< scalar_t > & HODLR_options()
Definition StrumpackOptions.hpp:1260
void set_proportional_mapping(ProportionalMapping pmap)
Definition StrumpackOptions.hpp:820
void set_components(int components)
Definition StrumpackOptions.hpp:418
Definition StrumpackOptions.hpp:44
GramSchmidtType
Definition StrumpackOptions.hpp:157
real_t default_rel_tol()
Definition StrumpackOptions.hpp:186
std::string get_description(MatchingJob job)
bool is_parallel(ReorderingStrategy method)
ProportionalMapping
Definition StrumpackOptions.hpp:51
MatchingJob
Definition StrumpackOptions.hpp:120
CompressionType
Definition StrumpackOptions.hpp:92
ReorderingStrategy
Definition StrumpackOptions.hpp:61
KrylovSolver
Definition StrumpackOptions.hpp:166
std::string get_name(ReorderingStrategy method)
MatchingJob get_matching(int job)
real_t default_abs_tol()
Definition StrumpackOptions.hpp:194