iMSTK
Interactive Medical Simulation Toolkit
imstkDirectLinearSolver.h
1 /*
2 ** This file is part of the Interactive Medical Simulation Toolkit (iMSTK)
3 ** iMSTK is distributed under the Apache License, Version 2.0.
4 ** See accompanying NOTICE for details.
5 */
6 
7 #pragma once
8 
9 #include "imstkLinearSolver.h"
10 #include "imstkMath.h"
11 
12 #ifdef WIN32
13 #pragma warning( push )
14 #pragma warning( disable : 4127 )
15 #endif
16 #include <Eigen/Sparse>
17 #include <Eigen/SparseLU>
18 #ifdef WIN32
19 #pragma warning( pop )
20 #endif
21 
22 namespace imstk
23 {
24 template<typename MatrixType> class DirectLinearSolver;
25 
32 template<>
33 class DirectLinearSolver<Matrixd>: public LinearSolver<Matrixd>
34 {
35 public:
36  DirectLinearSolver() { }
37  DirectLinearSolver(const Matrixd& A, const Vectord& b);
38  ~DirectLinearSolver() override { }
39 
43  void solve(Vectord& x) override;
44 
48  void solve(const Vectord& rhs, Vectord& x);
49 
53  void setSystem(std::shared_ptr<LinearSystemType> newSystem) override;
54 
58  void setSystem(std::shared_ptr<Matrixd> matrix);
59 
63  bool isIterative() const override
64  {
65  return false;
66  };
67 
68 private:
70  Eigen::LDLT<Matrixd> m_solver;
71 };
72 
77 template<>
78 class DirectLinearSolver<SparseMatrixd>: public LinearSolver<SparseMatrixd>
79 {
80 public:
84  DirectLinearSolver() = default;
85  ~DirectLinearSolver() override = default;
86 
90  DirectLinearSolver(const SparseMatrixd& matrix, const Vectord& b);
91 
95  void setSystem(std::shared_ptr<LinearSystemType> newSystem) override;
96 
100  void solve(Vectord& x) override;
101 
105  void solve(const Vectord& rhs, Vectord& x);
106 
107 private:
108  Eigen::SparseLU<SparseMatrixd, Eigen::COLAMDOrdering<MatrixType::StorageIndex>> m_solver;
109 };
110 } // namespace imstk
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > Matrixd
A dynamic size matrix of doubles.
Definition: imstkMath.h:97
Compound Geometry.
Dense direct solvers. Solves a dense system of equations using Cholesky decomposition.
Base class for linear solvers.
bool isIterative() const override
Returns true if the solver is iterative.