KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
eigen_dense_llt_solver.h
Go to the documentation of this file.
1 /* KRATOS _ _ ____ _
2 // | | (_)_ __ ___ __ _ _ __/ ___| ___ | |_ _____ _ __ ___
3 // | | | | '_ \ / _ \/ _` | '__\___ \ / _ \| \ \ / / _ \ '__/ __|
4 // | |___| | | | | __/ (_| | | ___) | (_) | |\ V / __/ | \__ |
5 // |_____|_|_| |_|\___|\__,_|_| |____/ \___/|_| \_/ \___|_| |___/ Application
6 //
7 // Author: Quirin Aumann
8 */
9 
10 #if !defined(KRATOS_EIGEN_DENSE_LLT_SOLVER_H_INCLUDED)
11 #define KRATOS_EIGEN_DENSE_LLT_SOLVER_H_INCLUDED
12 
13 // External includes
14 
15 // Project includes
16 #include "includes/define.h"
17 #include "linear_solvers_define.h"
18 
19 namespace Kratos {
20 
21 template <typename TScalar = double>
23 {
24 public:
25  using Scalar = TScalar;
28 
29 private:
30  Eigen::LLT<Matrix> m_solver;
31 
32 public:
33  static std::string Name()
34  {
35  return "complex_dense_llt";
36  }
37 
38  void Initialize(Parameters settings)
39  {
40  }
41 
42  bool Compute(Eigen::Map<Matrix> a)
43  {
44  m_solver.compute(a);
45 
46  const bool success = m_solver.info() == Eigen::Success;
47 
48  return success;
49  }
50 
51  bool Solve(Eigen::Ref<const Vector> b, Eigen::Ref<Vector> x) const
52  {
53  x = m_solver.solve(b);
54 
55  const bool success = m_solver.info() == Eigen::Success;
56 
57  return success;
58  }
59 
60  bool SolveMultiple(Eigen::Ref<const Matrix> b, Eigen::Ref<Matrix> x) const
61  {
62  x = m_solver.solve(b);
63 
64  const bool success = m_solver.info() == Eigen::Success;
65 
66  return success;
67  }
68 
69  void PrintInfo(std::ostream &rOStream) const
70  {
71  rOStream << "EigenDirectSolver <" << Name() << "> finished.";
72  }
73 
74  std::string GetSolverErrorMessages() const
75  {
76  return "No additional information";
77  }
78 };
79 
80 } // namespace Kratos
81 
82 #endif // defined(KRATOS_EIGEN_DENSE_LLT_SOLVER_H_INCLUDED)
Definition: eigen_dense_llt_solver.h:23
void Initialize(Parameters settings)
Definition: eigen_dense_llt_solver.h:38
static std::string Name()
Definition: eigen_dense_llt_solver.h:33
std::string GetSolverErrorMessages() const
Definition: eigen_dense_llt_solver.h:74
TScalar Scalar
Definition: eigen_dense_llt_solver.h:25
bool SolveMultiple(Eigen::Ref< const Matrix > b, Eigen::Ref< Matrix > x) const
Definition: eigen_dense_llt_solver.h:60
bool Compute(Eigen::Map< Matrix > a)
Definition: eigen_dense_llt_solver.h:42
bool Solve(Eigen::Ref< const Vector > b, Eigen::Ref< Vector > x) const
Definition: eigen_dense_llt_solver.h:51
Kratos::EigenDynamicMatrix< Scalar > Matrix
Definition: eigen_dense_llt_solver.h:26
Kratos::EigenDynamicVector< Scalar > Vector
Definition: eigen_dense_llt_solver.h:27
void PrintInfo(std::ostream &rOStream) const
Definition: eigen_dense_llt_solver.h:69
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Eigen::Matrix< _Scalar, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > EigenDynamicMatrix
Definition: linear_solvers_define.h:32
Eigen::Matrix< _Scalar, Eigen::Dynamic, 1 > EigenDynamicVector
Definition: linear_solvers_define.h:33
a
Definition: generate_stokes_twofluid_element.py:77
b
Definition: generate_total_lagrangian_mixed_volumetric_strain_element.py:31
x
Definition: sensitivityMatrix.py:49