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_pardiso_lu_solver.h
Go to the documentation of this file.
1 /* KRATOS _ _ ____ _
2 // | | (_)_ __ ___ __ _ _ __/ ___| ___ | |_ _____ _ __ ___
3 // | | | | '_ \ / _ \/ _` | '__\___ \ / _ \| \ \ / / _ \ '__/ __|
4 // | |___| | | | | __/ (_| | | ___) | (_) | |\ V / __/ | \__ |
5 // |_____|_|_| |_|\___|\__,_|_| |____/ \___/|_| \_/ \___|_| |___/ Application
6 //
7 // Author: Thomas Oberbichler
8 */
9 
10 #if !defined(KRATOS_EIGEN_PARDISO_LU_SOLVER_H_INCLUDED)
11 #define KRATOS_EIGEN_PARDISO_LU_SOLVER_H_INCLUDED
12 
13 // External includes
14 #include <Eigen/Sparse>
15 #include <Eigen/PardisoSupport>
16 
17 // Project includes
18 #include "includes/define.h"
19 #include "linear_solvers_define.h"
21 #include "spaces/ublas_space.h"
24 
25 namespace Kratos {
26 
27 template <typename TScalar = double>
29 {
30 public:
31  using Scalar = TScalar;
34 
35 private:
36  Eigen::PardisoLU<SparseMatrix> m_solver;
37 
38 public:
39  static std::string Name()
40  {
41  return "eigen_pardiso_lu";
42  }
43 
44  void Initialize(Parameters settings)
45  {
46  }
47 
48  bool Compute(Eigen::Map<const SparseMatrix> a)
49  {
50  m_solver.compute(a);
51 
52  const bool success = m_solver.info() == Eigen::Success;
53 
54  return success;
55  }
56 
57  bool Solve(Eigen::Ref<const Vector> b, Eigen::Ref<Vector> x) const
58  {
59  x = m_solver.solve(b);
60 
61  const bool success = m_solver.info() == Eigen::Success;
62 
63  return success;
64  }
65 
66  void PrintInfo(std::ostream &rOStream) const
67  {
68  rOStream << "EigenDirectSolver <" << Name() << "> finished.";
69  }
70 
71  std::string GetSolverErrorMessages() const
72  {
73  return "No additional information";
74  }
75 };
76 
77 } // namespace Kratos
78 
79 #endif // defined(KRATOS_EIGEN_PARDISO_LU_SOLVER_H_INCLUDED)
Definition: eigen_pardiso_lu_solver.h:29
Kratos::EigenDynamicVector< Scalar > Vector
Definition: eigen_pardiso_lu_solver.h:33
bool Solve(Eigen::Ref< const Vector > b, Eigen::Ref< Vector > x) const
Definition: eigen_pardiso_lu_solver.h:57
void Initialize(Parameters settings)
Definition: eigen_pardiso_lu_solver.h:44
Kratos::EigenSparseMatrix< Scalar > SparseMatrix
Definition: eigen_pardiso_lu_solver.h:32
std::string GetSolverErrorMessages() const
Definition: eigen_pardiso_lu_solver.h:71
bool Compute(Eigen::Map< const SparseMatrix > a)
Definition: eigen_pardiso_lu_solver.h:48
TScalar Scalar
Definition: eigen_pardiso_lu_solver.h:31
void PrintInfo(std::ostream &rOStream) const
Definition: eigen_pardiso_lu_solver.h:66
static std::string Name()
Definition: eigen_pardiso_lu_solver.h:39
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::SparseMatrix< _Scalar, Eigen::RowMajor, int > EigenSparseMatrix
Definition: linear_solvers_define.h:35
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