Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_Control.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Amesos2: Templated Direct Sparse Solver Package
4//
5// Copyright 2011 NTESS and the Amesos2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
17#ifndef AMESOS2_CONTROL_HPP
18#define AMESOS2_CONTROL_HPP
19
20#include <Teuchos_ParameterList.hpp>
21#include <Teuchos_RCP.hpp>
22
23namespace Amesos2 {
24
25
26struct Control {
28 Control()
29 : verbose_(0)
30 , debug_(0)
31 , useTranspose_(false)
32 , useIterRefine_(false)
33 , maxNumIterRefines_(2)
34 , verboseIterRefine_(false)
35 , addToDiag_("0.0")
36 , addZeroToDiag_(false)
37 , matrixProperty_(0)
38 , rcond_threshold_(1e-12)
39 , refactorize_(false)
40 , scaleMethod_(0)
41 , reindex_(false)
42 { }
43
45 ~Control() { };
46
47 void setControlParameters(
48 const Teuchos::RCP<Teuchos::ParameterList> & parameterList );
49
55 int verbose_;
56
58 int debug_;
59
60
62 bool useTranspose_;
63
65 bool useIterRefine_;
66
68 int maxNumIterRefines_;
69
71 bool verboseIterRefine_;
72
79 std::string addToDiag_;
80
81
87 bool addZeroToDiag_;
88
89
99 int matrixProperty_;
100
101
104 double rcond_threshold_; // if we refactorize, the factorization
105 // may suffer in numeric quality. We
106 // compute rcond = min (abs (diag (U))) /
107 // max (abs(diag (U))). If this ratio is
108 // <= rcond_threshold_, then the
109 // "refactorization" is scrapped, and we
110 // factor with full partial pivoting
111 // instead.
112
113
114 bool refactorize_; // if true, and if the Symbolic and Numeric
115 // objects have already been created, then
116 // attempt to "refactorize" (factor the matrix
117 // with no changes to the pivot order since the
118 // last call the klu_btf_factor).
119
120
121 int scaleMethod_; // Most methods (klu, UMFPACK, Mumps, ...) can
122 // scale the input matrix prior to
123 // factorization. This can improve pivoting,
124 // reduces fill-in, and leads to a better
125 // quality factorization. The options are:
126 // 0: no scaling
127 // 1: use the default method for the specific
128 // package
129 // 2: use the method's 1st alternative (if it
130 // has one)
131 // 3: use the method's 2nd alternative, and so
132 // on.
133 //
134 // Amesos2_Klu is, at present, the only code
135 // which implements this
136
137
145 bool reindex_;
146
147}; // end class Control
148
149
150} // end namespace Amesos2
151
152#endif // AMESOS2_CONTROL_HPP