15#ifndef dealii_trilinos_solver_h
16#define dealii_trilinos_solver_h
21#ifdef DEAL_II_WITH_TRILINOS
33# include <Epetra_LinearProblem.h>
34# include <Epetra_Operator.h>
37# ifdef DEAL_II_TRILINOS_WITH_BELOS
38# include <BelosBlockCGSolMgr.hpp>
39# include <BelosBlockGmresSolMgr.hpp>
40# include <BelosEpetraAdapter.hpp>
41# include <BelosIteration.hpp>
42# include <BelosMultiVec.hpp>
43# include <BelosOperator.hpp>
44# include <BelosSolverManager.hpp>
57 class PreconditionBase;
211 Epetra_MultiVector &x,
212 const Epetra_MultiVector &b,
227 Epetra_MultiVector &x,
228 const Epetra_MultiVector &b,
246 const ::Vector<double> &b,
263 const ::Vector<double> &b,
275 const ::LinearAlgebra::distributed::Vector<double> &b,
288 const ::LinearAlgebra::distributed::Vector<double> &b,
303 <<
"An error with error number " << arg1
304 <<
" occurred while calling a Trilinos function");
320 template <
typename Preconditioner>
322 do_solve(
const Preconditioner &preconditioner);
327 template <
typename Preconditioner>
577 const ::LinearAlgebra::distributed::Vector<double> &b);
594 const ::LinearAlgebra::distributed::Vector<double> &b)
const;
615 const ::Vector<double> &b);
626 const ::LinearAlgebra::distributed::Vector<double> &b);
639 <<
"An error with error number " << arg1
640 <<
" occurred while calling a Trilinos function");
674 std::unique_ptr<Amesos_BaseSolver>
solver;
684# ifdef DEAL_II_TRILINOS_WITH_BELOS
691 template <
typename VectorType>
717 struct AdditionalData
722 AdditionalData(
const SolverName solver_name = SolverName::cg,
723 const bool right_preconditioning =
false)
724 : solver_name(solver_name)
725 , right_preconditioning(right_preconditioning)
731 SolverName solver_name;
736 bool right_preconditioning;
743 const AdditionalData &additional_data,
744 const Teuchos::RCP<Teuchos::ParameterList> &belos_parameters);
749 template <
typename OperatorType,
typename PreconditionerType>
751 solve(
const OperatorType &a,
754 const PreconditionerType &p);
758 const AdditionalData additional_data;
759 const Teuchos::RCP<Teuchos::ParameterList> &belos_parameters;
769# ifdef DEAL_II_TRILINOS_WITH_BELOS
780 template <
typename VectorType>
782 class MultiVecWrapper
783 :
public Belos::MultiVec<typename VectorType::value_type>
789 using value_type =
typename VectorType::value_type;
795 this_type_is_missing_a_specialization()
803 MultiVecWrapper(VectorType &vector)
805 this->vectors.resize(1);
806 this->vectors[0].reset(
814 MultiVecWrapper(
const VectorType &vector)
816 this->vectors.resize(1);
817 this->vectors[0].reset(
825 virtual ~MultiVecWrapper() =
default;
830 virtual Belos::MultiVec<value_type> *
831 Clone(
const int numvecs)
const
833 auto new_multi_vec =
new MultiVecWrapper<VectorType>;
835 new_multi_vec->vectors.resize(numvecs);
837 for (
auto &vec : new_multi_vec->vectors)
839 vec = std::make_shared<VectorType>();
842 vec->reinit(*this->vectors[0]);
845 return new_multi_vec;
851 virtual Belos::MultiVec<value_type> *
862 virtual Belos::MultiVec<value_type> *
863 CloneCopy(
const std::vector<int> &index)
const
865 auto new_multi_vec =
new MultiVecWrapper<VectorType>;
867 new_multi_vec->vectors.resize(
index.size());
869 for (
unsigned int i = 0; i <
index.size(); ++i)
872 this->vectors.size(),
875 new_multi_vec->vectors[i] = std::make_shared<VectorType>();
878 *new_multi_vec->vectors[i] = *this->vectors[
index[i]];
881 return new_multi_vec;
889 virtual Belos::MultiVec<value_type> *
890 CloneViewNonConst(
const std::vector<int> &index)
892 auto new_multi_vec =
new MultiVecWrapper<VectorType>;
894 new_multi_vec->vectors.resize(
index.size());
896 for (
unsigned int i = 0; i <
index.size(); ++i)
899 this->vectors.size(),
902 new_multi_vec->vectors[i].reset(
903 this->vectors[index[i]].get(),
909 return new_multi_vec;
917 virtual const Belos::MultiVec<value_type> *
918 CloneView(
const std::vector<int> &index)
const
920 auto new_multi_vec =
new MultiVecWrapper<VectorType>;
922 new_multi_vec->vectors.resize(
index.size());
924 for (
unsigned int i = 0; i <
index.size(); ++i)
927 this->vectors.size(),
930 new_multi_vec->vectors[i].reset(
931 this->vectors[index[i]].get(),
937 return new_multi_vec;
944 GetGlobalLength()
const
948 for (
unsigned int i = 1; i < this->vectors.size(); ++i)
951 return this->vectors[0]->size();
958 GetNumberVecs()
const
960 return vectors.size();
967 MvTimesMatAddMv(
const value_type alpha,
968 const Belos::MultiVec<value_type> &A_,
969 const Teuchos::SerialDenseMatrix<int, value_type> &B,
970 const value_type beta)
972 const auto &A = try_to_get_underlying_vector(A_);
974 const unsigned int n_rows = B.numRows();
975 const unsigned int n_cols = B.numCols();
977 AssertThrow(n_rows ==
static_cast<unsigned int>(A.GetNumberVecs()),
979 AssertThrow(n_cols ==
static_cast<unsigned int>(this->GetNumberVecs()),
982 for (
unsigned int i = 0; i < n_cols; ++i)
983 (*this->vectors[i]) *= beta;
985 for (
unsigned int i = 0; i < n_cols; ++i)
986 for (
unsigned int j = 0; j < n_rows; ++j)
987 this->vectors[i]->add(alpha * B(j, i), *A.vectors[j]);
994 MvAddMv(
const value_type alpha,
995 const Belos::MultiVec<value_type> &A_,
996 const value_type beta,
997 const Belos::MultiVec<value_type> &B_)
999 const auto &A = try_to_get_underlying_vector(A_);
1000 const auto &B = try_to_get_underlying_vector(B_);
1002 AssertThrow(this->vectors.size() == A.vectors.size(),
1004 AssertThrow(this->vectors.size() == B.vectors.size(),
1007 for (
unsigned int i = 0; i < this->vectors.size(); ++i)
1009 this->vectors[i]->equ(alpha, *A.vectors[i]);
1010 this->vectors[i]->add(beta, *B.vectors[i]);
1018 MvScale(
const value_type alpha)
1020 for (
unsigned int i = 0; i < this->vectors.size(); ++i)
1021 (*this->vectors[i]) *= alpha;
1028 MvScale(
const std::vector<value_type> &alpha)
1039 MvTransMv(
const value_type alpha,
1040 const Belos::MultiVec<value_type> &A_,
1041 Teuchos::SerialDenseMatrix<int, value_type> &B)
const
1043 const auto &A = try_to_get_underlying_vector(A_);
1045 const unsigned int n_rows = B.numRows();
1046 const unsigned int n_cols = B.numCols();
1048 AssertThrow(n_rows ==
static_cast<unsigned int>(A.GetNumberVecs()),
1050 AssertThrow(n_cols ==
static_cast<unsigned int>(this->GetNumberVecs()),
1053 for (
unsigned int i = 0; i < n_rows; ++i)
1054 for (
unsigned int j = 0; j < n_cols; ++j)
1055 B(i, j) = alpha * ((*A.vectors[i]) * (*this->vectors[j]));
1063 MvDot(
const Belos::MultiVec<value_type> &A_,
1064 std::vector<value_type> &b)
const
1066 const auto &A = try_to_get_underlying_vector(A_);
1068 AssertThrow(this->vectors.size() == A.vectors.size(),
1072 for (
unsigned int i = 0; i < this->vectors.size(); ++i)
1073 b[i] = (*this->vectors[i]) * (*A.vectors[i]);
1081 std::vector<
typename Teuchos::ScalarTraits<value_type>::magnitudeType>
1083 Belos::NormType type = Belos::TwoNorm)
const
1088 for (
unsigned int i = 0; i < this->vectors.size(); ++i)
1089 normvec[i] = this->vectors[i]->
l2_norm();
1096 SetBlock(
const Belos::MultiVec<value_type> &A,
1097 const std::vector<int> &index)
1117 MvInit(
const value_type alpha)
1127 MvPrint(std::ostream &os)
const
1148 genericVector()
const
1158 static MultiVecWrapper<VectorType> &
1159 try_to_get_underlying_vector(Belos::MultiVec<value_type> &vec_in)
1161 auto vec =
dynamic_cast<MultiVecWrapper<VectorType> *
>(&vec_in);
1171 const static MultiVecWrapper<VectorType> &
1172 try_to_get_underlying_vector(
const Belos::MultiVec<value_type> &vec_in)
1174 auto vec =
dynamic_cast<const MultiVecWrapper<VectorType> *
>(&vec_in);
1182# ifdef HAVE_BELOS_TSQR
1184 factorExplicit(Belos::MultiVec<value_type> &,
1185 Teuchos::SerialDenseMatrix<int, value_type> &,
1193 Teuchos::SerialDenseMatrix<int, value_type> &,
1194 const typename Teuchos::ScalarTraits<value_type>::magnitudeType &)
1205 std::vector<std::shared_ptr<VectorType>> vectors;
1211 MultiVecWrapper() =
default;
1220 template <
typename OperatorType,
typename VectorType>
1222 class OperatorWrapper
1223 :
public Belos::Operator<typename VectorType::value_type>
1229 using value_type =
typename VectorType::value_type;
1235 this_type_is_missing_a_specialization()
1243 OperatorWrapper(
const OperatorType &op)
1250 virtual ~OperatorWrapper() =
default;
1256 Apply(
const Belos::MultiVec<value_type> &x,
1257 Belos::MultiVec<value_type> &y,
1258 Belos::ETrans trans = Belos::NOTRANS)
const
1263 op.vmult(MultiVecWrapper<VectorType>::try_to_get_underlying_vector(y)
1265 MultiVecWrapper<VectorType>::try_to_get_underlying_vector(x)
1273 HasApplyTranspose()
const
1283 const OperatorType &op;
1290 template <
typename VectorType>
1292 SolverBelos<VectorType>::SolverBelos(
1294 const AdditionalData &additional_data,
1295 const Teuchos::RCP<Teuchos::ParameterList> &belos_parameters)
1296 : solver_control(solver_control)
1297 , additional_data(additional_data)
1298 , belos_parameters(belos_parameters)
1303 template <
typename VectorType>
1305 template <
typename OperatorType,
typename PreconditionerType>
1306 void SolverBelos<VectorType>::solve(
const OperatorType &A_dealii,
1307 VectorType &x_dealii,
1308 const VectorType &b_dealii,
1309 const PreconditionerType &P_dealii)
1311 using value_type =
typename VectorType::value_type;
1313 using MV = Belos::MultiVec<value_type>;
1314 using OP = Belos::Operator<value_type>;
1316 Teuchos::RCP<OP> A = Teuchos::rcp(
1317 new internal::OperatorWrapper<OperatorType, VectorType>(A_dealii));
1318 Teuchos::RCP<OP> P = Teuchos::rcp(
1319 new internal::OperatorWrapper<PreconditionerType, VectorType>(P_dealii));
1320 Teuchos::RCP<MV> X =
1321 Teuchos::rcp(
new internal::MultiVecWrapper<VectorType>(x_dealii));
1322 Teuchos::RCP<MV> B =
1323 Teuchos::rcp(
new internal::MultiVecWrapper<VectorType>(b_dealii));
1325 Teuchos::RCP<Belos::LinearProblem<value_type, MV, OP>> problem =
1326 Teuchos::rcp(
new Belos::LinearProblem<value_type, MV, OP>(A, X, B));
1328 if (additional_data.right_preconditioning ==
false)
1329 problem->setLeftPrec(P);
1331 problem->setRightPrec(P);
1333 const bool problem_set = problem->setProblem();
1338 r.reinit(x_dealii,
true);
1339 A_dealii.vmult(r, x_dealii);
1340 r.sadd(-1., 1., b_dealii);
1341 const auto norm_0 = r.l2_norm();
1346 double relative_tolerance_to_be_achieved =
1347 solver_control.tolerance() / norm_0;
1348 const unsigned int max_steps = solver_control.max_steps();
1350 if (
const auto *reduction_control =
1352 relative_tolerance_to_be_achieved =
1353 std::max(relative_tolerance_to_be_achieved,
1354 reduction_control->reduction());
1356 Teuchos::RCP<Teuchos::ParameterList> belos_parameters_copy(
1357 Teuchos::rcp(
new Teuchos::ParameterList(*belos_parameters)));
1359 belos_parameters_copy->set(
"Convergence Tolerance",
1360 relative_tolerance_to_be_achieved);
1361 belos_parameters_copy->set(
"Maximum Iterations",
1362 static_cast<int>(max_steps));
1364 Teuchos::RCP<Belos::SolverManager<value_type, MV, OP>> solver;
1366 if (additional_data.solver_name == SolverName::cg)
1367 solver = Teuchos::rcp(
1368 new Belos::BlockCGSolMgr<value_type, MV, OP>(problem,
1369 belos_parameters_copy));
1370 else if (additional_data.solver_name == SolverName::gmres)
1371 solver = Teuchos::rcp(
1372 new Belos::BlockGmresSolMgr<value_type, MV, OP>(problem,
1373 belos_parameters_copy));
1377 const auto flag = solver->solve();
1379 solver_control.check(solver->getNumIters(), solver->achievedTol() * norm_0);
1381 AssertThrow(flag == Belos::ReturnType::Converged ||
1384 (solver_control.last_step() == max_steps)),
1386 solver_control.last_value()));
@ iterate
Continue iteration.
virtual ~SolverBase()=default
SolverControl & control() const
void solve(const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner)
const AdditionalData additional_data
std::unique_ptr< AztecOO_StatusTest > status_test
void set_preconditioner(AztecOO &solver, const Preconditioner &preconditioner)
SolverControl & solver_control
void do_solve(const Preconditioner &preconditioner)
void solve(Epetra_Operator &A, ::LinearAlgebra::distributed::Vector< double > &x, const ::LinearAlgebra::distributed::Vector< double > &b, const PreconditionBase &preconditioner)
std::unique_ptr< Epetra_LinearProblem > linear_problem
enum TrilinosWrappers::SolverBase::SolverName solver_name
void solve(const SparseMatrix &A, ::LinearAlgebra::distributed::Vector< double > &x, const ::LinearAlgebra::distributed::Vector< double > &b, const PreconditionBase &preconditioner)
void solve(const SparseMatrix &A, ::LinearAlgebra::distributed::Vector< double > &x, const ::LinearAlgebra::distributed::Vector< double > &b)
virtual ~SolverDirect()=default
void initialize(const SparseMatrix &A)
void vmult(MPI::Vector &x, const MPI::Vector &b) const
std::unique_ptr< Amesos_BaseSolver > solver
SolverControl & solver_control
void solve(MPI::Vector &x, const MPI::Vector &b)
std::unique_ptr< Epetra_LinearProblem > linear_problem
AdditionalData additional_data
SolverControl solver_control_own
SolverControl & control() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
#define DeclException1(Exception1, type1, outsequence)
#define AssertThrow(cond, exc)
#define DEAL_II_NOT_IMPLEMENTED()
std::vector< index_type > data
std::vector< value_type > l2_norm(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
Tpetra::Vector< Number, LO, GO, NodeType< MemorySpace > > VectorType
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
const unsigned int gmres_restart_parameter
const bool output_solver_details
bool output_solver_details