HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
liealgebra.hpp
Go to the documentation of this file.
1#ifndef LIE_ALGEBRA
2#define LIE_ALGEBRA
3
4#include <Eigen/Dense>
5#include <vector>
6
7namespace HArDCore3D
8{
9
17 {
18 public:
19
20 typedef Eigen::MatrixXcd LieAlgValue;
21 typedef std::function<double(LieAlgValue &, LieAlgValue &)> LieAlgProduct;
22
24 LieAlgebra();
25 LieAlgebra(std::vector<LieAlgValue> & basis, LieAlgProduct & product);
26
28 inline const size_t dimension() const
29 {
30 return m_basis.size();
31 }
32
34 inline const std::vector<LieAlgValue> & basis() const
35 {
36 return m_basis;
37 }
38
40 inline const Eigen::MatrixXd & massMatrix() const
41 {
42 return m_mass_matrix;
43 }
44
46 inline const LieAlgValue lieBracket(const LieAlgValue & A, const LieAlgValue & B) const
47 {
48 return A*B - B*A;
49 }
51 inline const std::vector<Eigen::MatrixXd> & structureConst() const
52 {
53 return m_strucConst;
54 }
55
56 private:
57
59 Eigen::MatrixXd _compute_mass_matrix();
61 std::vector<Eigen::MatrixXd> _compute_structure_constants();
62
63 std::vector<LieAlgValue> m_basis;
64 LieAlgProduct m_product;
65 Eigen::MatrixXd m_mass_matrix;
66 std::vector<Eigen::MatrixXd> m_strucConst;
67 };
68
69}
70
71#endif
Lie algebra class: mass matrix, structure constants and Lie bracket.
Definition liealgebra.hpp:17
@ Matrix
Definition basis.hpp:67
std::function< double(LieAlgValue &, LieAlgValue &)> LieAlgProduct
Definition liealgebra.hpp:21
LieAlgebra()
Constructors.
Definition liealgebra.cpp:10
Eigen::MatrixXcd LieAlgValue
Definition liealgebra.hpp:20
const std::vector< Eigen::MatrixXd > & structureConst() const
Computes the structure constants of the Lie algebra.
Definition liealgebra.hpp:51
const std::vector< LieAlgValue > & basis() const
Returns the basis of the Lie algebra.
Definition liealgebra.hpp:34
const Eigen::MatrixXd & massMatrix() const
Returns the Gram matrix of the Lie algebra.
Definition liealgebra.hpp:40
const size_t dimension() const
Assuming orthonormal basis.
Definition liealgebra.hpp:28
const LieAlgValue lieBracket(const LieAlgValue &A, const LieAlgValue &B) const
Computes the Lie bracket of two Lie algebra elements.
Definition liealgebra.hpp:46
Definition ddr-magnetostatics.hpp:41