HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Loading...
Searching...
No Matches
xgrad-test.hpp
Go to the documentation of this file.
1#ifndef DDRSPACES_TEST_HPP
2#define DDRSPACES_TEST_HPP
3
4#include <boost/math/constants/constants.hpp>
5
6namespace HArDCore2D
7{
8
9 static const double PI = boost::math::constants::pi<double>();
10 using std::sin;
11
12 //------------------------------------------------------------------------------
13
14 static std::function<double(const Eigen::Vector2d&)>
15 trigonometric_scalar = [](const Eigen::Vector2d & x) -> double {
16 return sin(PI * x(0)) * sin(PI * x(1));
17 };
18
19 static std::function<Eigen::Vector2d(const Eigen::Vector2d&)>
20 grad_trigonometric_scalar = [](const Eigen::Vector2d & x) -> Eigen::Vector2d {
21 return PI * Eigen::Vector2d(
22 cos(PI * x(0)) * sin(PI * x(1)),
23 sin(PI * x(0)) * cos(PI * x(1))
24 );
25 };
26 //------------------------------------------------------------------------------
27
28 static std::function<double(const Eigen::Vector2d&)>
29 constant_scalar = [](const Eigen::Vector2d & x) -> double {
30 return 1.;
31 };
32
33 static std::function<Eigen::Vector2d(const Eigen::Vector2d&)>
34 grad_constant_scalar = [](const Eigen::Vector2d & x) -> Eigen::Vector2d {
35 return Eigen::Vector2d::Zero();
36 };
37
38 //------------------------------------------------------------------------------
39
40 static std::function<double(const Eigen::Vector2d&)>
41 linear_scalar = [](const Eigen::Vector2d & x) -> double {
42 return 1. + x(0) + 2. * x(1);
43 };
44
45 static std::function<Eigen::Vector2d(const Eigen::Vector2d&)>
46 grad_linear_scalar = [](const Eigen::Vector2d & x) -> Eigen::Vector2d {
47 return Eigen::Vector2d(1., 2.);
48 };
49
50 //------------------------------------------------------------------------------
51
52 static std::function<double(const Eigen::Vector2d&)>
53 quadratic_scalar = [](const Eigen::Vector2d & x) -> double {
54 return linear_scalar(x) + std::pow(x(0), 2) + 2. * std::pow(x(1), 2);
55 };
56
57 static std::function<Eigen::Vector2d(const Eigen::Vector2d&)>
58 grad_quadratic_scalar = [](const Eigen::Vector2d & x) -> Eigen::Vector2d {
59 return grad_linear_scalar(x) + Eigen::Vector2d(2. * x(0), 4. * x(1));
60 };
61
62 //------------------------------------------------------------------------------
63
64 template<typename T>
65 double squared_l2_error(
66 const std::function<T(const Eigen::Vector2d &)> & f,
67 const Eigen::VectorXd & fX,
68 const boost::multi_array<T, 2> & fX_basis_quad,
69 const QuadratureRule & quad_X
70 )
71 {
72 // Check that the dimensions are consistent
73 assert(
74 fX_basis_quad.shape()[0] == (size_t)fX.size() &&
75 fX_basis_quad.shape()[1] == quad_X.size()
76 );
77
78 double err = 0.;
79
80 for (size_t iqn = 0; iqn < quad_X.size(); iqn++) {
81 T f_iqn = f(quad_X[iqn].vector());
82
83 T fX_iqn = fX(0) * fX_basis_quad[0][iqn];
84 for (size_t i = 1; i < fX_basis_quad.shape()[0]; i++) {
85 fX_iqn += fX(i) * fX_basis_quad[i][iqn];
86 } // for i
87
88 T diff_iqn = f_iqn - fX_iqn;
89
90 err += quad_X[iqn].w * scalar_product(diff_iqn, diff_iqn);
91 } // for iqn
92
93 return err;
94 }
95
96
97} // end of namespace HArDCore2D
98
99#endif
Compute max and min eigenvalues of all matrices for i
Definition compute_eigs.m:5
double scalar_product(const double &x, const double &y)
Scalar product between two reals.
Definition basis.cpp:163
static const double PI
Definition ddr-klplate.hpp:187
std::vector< QuadratureNode > QuadratureRule
Definition quadraturerule.hpp:55
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
Definition ddr-klplate.hpp:27
static std::function< double(const Eigen::Vector2d &)> constant_scalar
Definition xgrad-test.hpp:29
static std::function< double(const Eigen::Vector2d &)> quadratic_scalar
Definition xgrad-test.hpp:53
static std::function< Eigen::Vector2d(const Eigen::Vector2d &)> grad_quadratic_scalar
Definition xgrad-test.hpp:58
static std::function< Eigen::Vector2d(const Eigen::Vector2d &)> grad_trigonometric_scalar
Definition xgrad-test.hpp:20
static std::function< double(const Eigen::Vector2d &)> trigonometric_scalar
Definition xgrad-test.hpp:15
static std::function< double(const Eigen::Vector2d &)> linear_scalar
Definition xgrad-test.hpp:41
static std::function< Eigen::Vector2d(const Eigen::Vector2d &)> grad_linear_scalar
Definition xgrad-test.hpp:46
static std::function< Eigen::Vector2d(const Eigen::Vector2d &)> grad_constant_scalar
Definition xgrad-test.hpp:34
double squared_l2_error(const std::function< T(const Eigen::Vector2d &)> &f, const Eigen::VectorXd &fX, const boost::multi_array< T, 2 > &fX_basis_quad, const QuadratureRule &quad_X)
Definition excurl-test.hpp:140