HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
bgg-klplate.hpp
Go to the documentation of this file.
1// Implementation of the DDR-BGG scheme for the Kirchhoff-Love plate problem.
2//
3// Author: Arax Leroy (arax.leroy@umontpellier.fr)
4//
5
6
7#ifndef BGG_KLPLATE_HPP
8#define BGG_KLPLATE_HPP
9
10#include <iostream>
11
12#include <boost/math/constants/constants.hpp>
13
14#include <Eigen/Sparse>
15#include <unsupported/Eigen/SparseExtra>
16
17#include <mesh.hpp>
18#include <BoundaryConditions/BoundaryConditions.hpp> // To re-order the boundary edges and vertices
20
21#include <vsxgrad.hpp>
22#include <xhess.hpp>
23
29namespace HArDCore2D
30{
31
38 struct KLNorms
39 {
41 KLNorms(double norm_displacement):
42 displacement(norm_displacement),
43 energy(std::sqrt( std::pow(norm_displacement, 2) ) )
44 {
45 // do nothing
46 };
47 double displacement;
48 double energy;
49 };
50
53 {
54 typedef Eigen::SparseMatrix<double> SystemMatrixType;
55 typedef std::function<double(const Eigen::Vector2d &)> ForcingTermType;
56 typedef std::function<double(const Eigen::Vector2d &)> SolutionDisplacementType;
57 typedef std::function<Eigen::Vector2d(const Eigen::Vector2d&)> GradientDisplacementType;
58
59
62 const DDRCore & ddrcore,
63 const DDRCore & ddrcore_plus,
64 const SerendipityProblem & sp,
65 const SerendipityProblem & sp_plus,
66 const BoundaryConditions & BC_u,
67 bool use_threads,
68 std::ostream & output = std::cout
69 );
70
73 const ForcingTermType & f,
74 const SolutionDisplacementType & u,
75 const GradientDisplacementType & grad_u
76 );
77
79 inline size_t dimensionSpace() const
80 {
81 return m_xhess.dimension();
82 }
83
85 inline size_t nb_bdryDOFs() const
86 {
87 return m_BC_u.n_dir_edges() * m_xhess.numLocalDofsEdge()
88 + m_BC_u.n_dir_vertices() * m_xhess.numLocalDofsVertex();
89 }
90
92 inline size_t sizeSystem() const
93 {
94 return dimensionSpace() - nb_bdryDOFs();
95 }
96
98 inline const std::vector<std::pair<size_t,size_t>> & locUKN() const {
99 return m_locUKN;
100 }
101
103 std::vector<size_t> globalDOFIndices(const Cell &T) const
104 {
105 return m_xhess.globalDOFIndices(T);
106 }
107
109 inline const VSXGrad & vsxGrad() const
110 {
111 return m_vsxgrad;
112 }
113
115 inline const XHess & xHess() const
116 {
117 return m_xhess;
118 }
119
121 inline const SystemMatrixType & systemMatrix() const {
122 return m_A;
123 }
124
127 return m_A;
128 }
129
131 inline const Eigen::VectorXd & systemVector() const {
132 return m_b;
133 }
134
136 inline Eigen::VectorXd & systemVector() {
137 return m_b;
138 }
139
141 inline const SystemMatrixType & bdryMatrix() const {
142 return m_bdryMatrix;
143 }
144
146 inline const Eigen::VectorXd & bdryValues() const {
147 return m_bdryValues;
148 }
149
151 inline const double & stabilizationParameter() const {
152 return m_stab_par;
153 }
154
156 inline double & stabilizationParameter() {
157 return m_stab_par;
158 }
159
162 const Eigen::VectorXd & v
163 ) const;
164
166 template<typename outValue, typename Fct>
167 std::function<outValue(const Eigen::Vector2d &)> contractPara(const Fct &F) const
168 {
169 std::function<outValue(const Eigen::Vector2d &)> f = [this, &F](const Eigen::Vector2d &x)->outValue { return F(x);};
170 return f;
171 }
172
173 private:
175 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
176 _compute_local_contribution(
177 size_t iT,
178 const ForcingTermType & f
179 );
180
181
183
185 void _assemble_local_contribution(
186 size_t iT,
187 const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
188 std::list<Eigen::Triplet<double> > & A1,
189 Eigen::VectorXd & b1,
190 std::list<Eigen::Triplet<double> > & A2
191 );
192
193 const DDRCore & m_ddrcore;
194 bool m_use_threads;
195 std::ostream & m_output;
196 VSXGrad m_vsxgrad;
197 XHess m_xhess;
198 BoundaryConditions m_BC_u;
199 SystemMatrixType m_bdryMatrix; // Matrix with columns corresponding to Dirichlet DOFs; multiplied by the boundary values, yields the modification to the system RHS for non-homogeneous BC
200 Eigen::VectorXd m_bdryValues; // Boundary values
201 SystemMatrixType m_A; // Matrix and RHS for system (after removing Dirichlet DOFs)
202 Eigen::VectorXd m_b;
203 double m_stab_par;
204 std::vector<std::pair<size_t,size_t>> m_locUKN; // Indicates the location, among the DOFs, of the unknowns after the Dirichlet DOFs are removed (the unknowns are between m_locUCK[i].first and m_locUCK[i].first+m_locUKN[i].second for all i)
205 Eigen::VectorXi m_DOFtoUKN; // Maps a dof i to its unknown in the system without BC, or to -1 if the DOF is a Dirichlet DOF
206 };
207
208 //------------------------------------------------------------------------------
209 // Exact solutions
210 //------------------------------------------------------------------------------
211
212 static const double PI = boost::math::constants::pi<double>();
213 using std::sin;
214 using std::cos;
215 using std::exp;
216 using std::pow;
217
218 //------------------------------------------------------------------------------
219 // Constant
221 constant_u = [](const VectorRd & x) -> double {
222 return 1.;
223 };
224
226 constant_grad_u = [](const VectorRd & /*x*/) -> VectorRd{
227 return VectorRd::Zero();
228 };
229
231 constant_f = [](const VectorRd & x) -> double {
232 return 0.;
233 };
234
235 //------------------------------------------------------------------------------
236 // Polynomial
237
239 polynomial_u = [](const VectorRd & x) -> double {
240 return pow(x(0), 5) + pow(x(1), 5);
241 };
242
245 VectorRd grad;
246 grad(0) = 5.0 * pow(x(0), 4);
247 grad(1) = 5.0 * pow(x(1), 4);
248 return grad;
249 };
250
252 polynomial_f = [](const VectorRd & x) -> double {
253 return 120. * (x(0) + x(1));
254 };
255
256} // end of namespace HArDCore2D
257
258#endif
The BoundaryConditions class provides definition of boundary conditions.
Definition BoundaryConditions.hpp:45
const size_t n_dir_vertices() const
Returns the number of Dirichlet vertices.
Definition BoundaryConditions.hpp:75
const size_t n_dir_edges() const
Returns the number of Dirichlet edges.
Definition BoundaryConditions.hpp:70
Construct all polynomial spaces for the DDR sequence.
Definition ddrcore.hpp:63
Construct all polynomial spaces for the DDR sequence.
Definition serendipity_problem.hpp:20
Vector version of sXgrad, the arbitrary order space with nodal primal unknowns.
Definition vsxgrad.hpp:32
Discrete H2 space: local operators, L2 product and global interpolator.
Definition xhess.hpp:21
Create grid points x
Definition generate_cartesian_mesh.m:22
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition bgg-klplate.hpp:136
std::function< Eigen::Vector2d(const Eigen::Vector2d &)> GradientDisplacementType
Definition bgg-klplate.hpp:57
const Eigen::VectorXd & bdryValues() const
Returns the boundary values.
Definition bgg-klplate.hpp:146
double & stabilizationParameter()
Returns the stabilization parameter.
Definition bgg-klplate.hpp:156
const std::vector< std::pair< size_t, size_t > > & locUKN() const
Returns the location of the unknowns among the DOFs.
Definition bgg-klplate.hpp:98
KLNorms computeNorms(const Eigen::VectorXd &v) const
Compute the discrete norms: displacement,and energy.
Definition bgg-klplate.cpp:467
std::vector< size_t > globalDOFIndices(const Cell &T) const
Create the vector of DOF indices for cell T, which combines the DOFs for the spaces VSXgrad and Xhess...
Definition bgg-klplate.hpp:103
const SystemMatrixType & bdryMatrix() const
Returns the Matrix for BC.
Definition bgg-klplate.hpp:141
size_t nb_bdryDOFs() const
Returns the nb of DOFs for BC.
Definition bgg-klplate.hpp:85
static KirchhoffLove::SolutionDisplacementType constant_u
Definition bgg-klplate.hpp:221
std::function< outValue(const Eigen::Vector2d &)> contractPara(const Fct &F) const
Takes a function dependent on KLParameter and a position x, and returns a function depending only on ...
Definition bgg-klplate.hpp:167
static KirchhoffLove::GradientDisplacementType polynomial_grad_u
Definition bgg-klplate.hpp:244
const XHess & xHess() const
Returns the space XHess.
Definition bgg-klplate.hpp:115
void assembleLinearSystem(const ForcingTermType &f, const SolutionDisplacementType &u, const GradientDisplacementType &grad_u)
Assemble the global system
Definition bgg-klplate.cpp:292
Eigen::SparseMatrix< double > SystemMatrixType
Definition bgg-klplate.hpp:54
static const double PI
Definition bgg-klplate.hpp:212
std::function< double(const Eigen::Vector2d &)> ForcingTermType
Definition bgg-klplate.hpp:55
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition bgg-klplate.hpp:121
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition bgg-klplate.hpp:131
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition bgg-klplate.hpp:126
static KirchhoffLove::ForcingTermType constant_f
Definition bgg-klplate.hpp:231
static KirchhoffLove::ForcingTermType polynomial_f
Definition bgg-klplate.hpp:252
std::function< double(const Eigen::Vector2d &)> SolutionDisplacementType
Definition bgg-klplate.hpp:56
double displacement
Norm of displacement.
Definition bgg-klplate.hpp:47
size_t sizeSystem() const
Returns the size of the system without BC.
Definition bgg-klplate.hpp:92
KLNorms(double norm_displacement)
Constructor.
Definition bgg-klplate.hpp:41
const VSXGrad & vsxGrad() const
Returns the space VSXGrad.
Definition bgg-klplate.hpp:109
static KirchhoffLove::GradientDisplacementType constant_grad_u
Definition bgg-klplate.hpp:226
double energy
Total energy.
Definition bgg-klplate.hpp:48
const double & stabilizationParameter() const
Returns the stabilization parameter.
Definition bgg-klplate.hpp:151
static KirchhoffLove::SolutionDisplacementType polynomial_u
Definition bgg-klplate.hpp:239
size_t dimensionSpace() const
Returns the dimension of the rotation + displacement space (with BC)
Definition bgg-klplate.hpp:79
Eigen::Vector2d VectorRd
Definition basis.hpp:55
std::vector< size_t > globalDOFIndices(const Cell &T) const
Definition globaldofspace.cpp:98
size_t numLocalDofsVertex() const
Returns the number of local vertex DOFs.
Definition localdofspace.hpp:39
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition localdofspace.hpp:61
size_t numLocalDofsEdge() const
Returns the number of local edge DOFs.
Definition localdofspace.hpp:45
bool use_threads
Definition HHO_DiffAdvecReac.hpp:45
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
Definition bgg-klplate.hpp:30
static auto v
Definition ddrcore-test.hpp:32
Structure to store component norms (displacement total energy)
Definition bgg-klplate.hpp:39
Assemble a RM problem.
Definition bgg-klplate.hpp:53