HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
xrot.hpp
Go to the documentation of this file.
1#ifndef XROT_HPP
2#define XROT_HPP
3
4#include <globaldofspace.hpp>
5#include <ddrcore.hpp>
6#include <integralweight.hpp>
7#include <xrotrot.hpp>
8
9namespace HArDCore2D
10{
17 class XRot : public GlobalDOFSpace
18 {
19 public:
20 typedef std::function<double(const Eigen::Vector2d &)> FunctionType;
21
24 {
26 const Eigen::MatrixXd & _rotor,
27 const Eigen::MatrixXd & _rotor_rhs,
28 const Eigen::MatrixXd & _potential
29 )
30 : rotor(_rotor),
31 rotor_rhs(_rotor_rhs),
32 potential(_potential)
33 {
34 // Do nothing
35 }
36
37 Eigen::MatrixXd rotor;
38 Eigen::MatrixXd rotor_rhs;
39 Eigen::MatrixXd potential;
40 };
41
43 XRot(const DDRCore & ddr_core, bool use_threads = true, std::ostream & output = std::cout);
44
46 const Mesh & mesh() const
47 {
48 return m_ddr_core.mesh();
49 }
50
52 const bool useThreads() const
53 {
54 return m_use_threads;
55 }
56
58 const size_t & degree() const
59 {
60 return m_ddr_core.degree();
61 }
62
64 Eigen::VectorXd interpolate(
65 const FunctionType & q,
66 const int deg_quad = -1
67 ) const;
68
70 inline const Eigen::MatrixXd & edgePotential(size_t iE) const
71 {
72 return *m_edge_potentials[iE];
73 }
74
76 inline const Eigen::MatrixXd & edgePotential(const Edge& E) const
77 {
78 return *m_edge_potentials[E.global_index()];
79 }
80
82 inline const LocalOperators & cellOperators(size_t iT) const
83 {
84 return *m_cell_operators[iT];
85 }
86
88 inline const LocalOperators & cellOperators(const Cell & T) const
89 {
90 return *m_cell_operators[T.global_index()];
91 }
92
94 inline const DDRCore::CellBases & cellBases(size_t iT) const
95 {
96 return m_ddr_core.cellBases(iT);
97 }
98
100 inline const DDRCore::CellBases & cellBases(const Cell & T) const
101 {
102 return m_ddr_core.cellBases(T.global_index());
103 }
104
106 inline const DDRCore::EdgeBases & edgeBases(size_t iE) const
107 {
108 return m_ddr_core.edgeBases(iE);
109 }
110
112 inline const DDRCore::EdgeBases & edgeBases(const Edge & E) const
113 {
114 return m_ddr_core.edgeBases(E.global_index());
115 }
116
118 Eigen::MatrixXd computeRotorL2Product(
119 const size_t iT,
120 const double & penalty_factor_cell = 1.,
121 const double & penalty_factor_edge = 1.
122
123 ) const;
124
126 double computeRotorL2Norm(const Eigen::VectorXd & v) const;
127
128 private:
129 Eigen::MatrixXd _compute_edge_potential(size_t iE);
130 LocalOperators _compute_cell_operators(size_t iT);
131
132 const DDRCore & m_ddr_core;
133 bool m_use_threads;
134 std::ostream & m_output;
135
136 // Containers for local operators
137 std::vector<std::unique_ptr<Eigen::MatrixXd> > m_edge_potentials;
138 std::vector<std::unique_ptr<LocalOperators> > m_cell_operators;
139 };
140
141} // end of namespace HArDCore2D
142
143#endif
Construct all polynomial spaces for the DDR sequence.
Definition ddrcore.hpp:63
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition globaldofspace.hpp:16
Discrete H1 space: local operators, L2 product and global interpolator.
Definition xrot.hpp:18
Definition Mesh2D.hpp:26
const Mesh & mesh() const
Return a const reference to the mesh.
Definition ddrcore.hpp:115
const bool useThreads() const
Return true if we use thread-based parallelism.
Definition xrot.hpp:52
Eigen::MatrixXd computeRotorL2Product(const size_t iT, const double &penalty_factor_cell=1., const double &penalty_factor_edge=1.) const
Compute rotor L2-product.
Definition xrot.cpp:287
const size_t & degree() const
Return the polynomial degree.
Definition ddrcore.hpp:121
Eigen::MatrixXd rotor_rhs
Definition xrot.hpp:38
const DDRCore::EdgeBases & edgeBases(size_t iE) const
Return edge bases for the edge of index iE.
Definition xrot.hpp:106
const size_t & degree() const
Return the polynomial degree.
Definition xrot.hpp:58
double computeRotorL2Norm(const Eigen::VectorXd &v) const
Compute rotor L2-norm.
Definition xrot.cpp:348
Eigen::MatrixXd rotor
Definition xrot.hpp:37
const DDRCore::CellBases & cellBases(const Cell &T) const
Return cell bases for cell T.
Definition xrot.hpp:100
std::function< double(const Eigen::Vector2d &)> FunctionType
Definition xrot.hpp:20
Eigen::MatrixXd potential
Definition xrot.hpp:39
const DDRCore::CellBases & cellBases(size_t iT) const
Return cell bases for the cell of index iT.
Definition xrot.hpp:94
const CellBases & cellBases(size_t iT) const
Return cell bases for element iT.
Definition ddrcore.hpp:127
const Eigen::MatrixXd & edgePotential(const Edge &E) const
Return edge potential for the edge E.
Definition xrot.hpp:76
Eigen::VectorXd interpolate(const FunctionType &q, const int deg_quad=-1) const
Interpolator of a continuous function.
Definition xrot.cpp:61
const Mesh & mesh() const
Return the mesh.
Definition xrot.hpp:46
const DDRCore::EdgeBases & edgeBases(const Edge &E) const
Return edge bases for edge E.
Definition xrot.hpp:112
const Eigen::MatrixXd & edgePotential(size_t iE) const
Return edge potential for the edge of index iE.
Definition xrot.hpp:70
const EdgeBases & edgeBases(size_t iE) const
Return edge bases for edge iE.
Definition ddrcore.hpp:135
const LocalOperators & cellOperators(const Cell &T) const
Return cell operators for cell T.
Definition xrot.hpp:88
LocalOperators(const Eigen::MatrixXd &_rotor, const Eigen::MatrixXd &_rotor_rhs, const Eigen::MatrixXd &_potential)
Definition xrot.hpp:25
const LocalOperators & cellOperators(size_t iT) const
Return cell operators for the cell of index iT.
Definition xrot.hpp:82
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
Definition ddr-klplate.hpp:27
static auto v
Definition ddrcore-test.hpp:32
static auto q
Definition ddrcore-test.hpp:14
Structure to store element bases.
Definition ddrcore.hpp:81
Structure to store edge bases.
Definition ddrcore.hpp:102
A structure to store local operators (vector rotor and potential)
Definition xrot.hpp:24