Routines for the evalutation of pointwise errors. More...
Routines for the evalutation of pointwise errors.
| double Bembel::calculateFirstCoefficient | ( | Eigen::VectorXd | cs, |
| unsigned int | deg, | ||
| Eigen::MatrixXd | ps_l, | ||
| Eigen::MatrixXd | ps_f, | ||
| Eigen::MatrixXd | ps_b | ||
| ) |
Calculates the first coefficient such that the mean of the kernel vanishes.
| cs | The already calculated other coefficients |
| deg | The degree |
| ps_l | Gauss quadrature points on the left side of the cube |
| ps_f | Gauss quadrature points on the front side of the cube |
| ps_b | Gauss quadrature points on the bottom side of the cube |
Definition at line 346 of file Coefficients.hpp.
| std::string Bembel::double_to_string | ( | double | d, |
| const int | precision | ||
| ) |
Convert a double to a string with given precision.
| d | double to be converted. |
| precision | Precision of the conversion. |
Definition at line 253 of file GeometryIGS.hpp.
| Eigen::Vector3d Bembel::evaluate_dsolid_sphericals | ( | Eigen::Vector3d | x, |
| Eigen::VectorXd | cs, | ||
| unsigned int | deg | ||
| ) |
Evaluates the series \( \sum_{n = 0}^{\rm deg} \sum_{m = -n}^n |x|^{n-3} c_m^n (\nabla Y_m^n)(\frac{x}{|x|}) \) for real coefficients.
| x | The point of evaluation |
| cs | The coefficients stored in the order \( [(0, 0), (1, -1), (1, 0), (1, 1) (2, -2), (2, -1), ... (n, n)] \) |
| deg | The degree |
Definition at line 302 of file Sphericals.hpp.
| Eigen::Vector3d Bembel::evaluate_dsphericals | ( | Eigen::Vector3d | x, |
| Eigen::VectorXd | cs, | ||
| unsigned int | deg | ||
| ) |
Evaluates the series \( \sum_{n = 0}^{\rm deg} \sum_{m = -n}^n c_m^n \nabla Y_m^n(x) \) for real coefficients.
| x | The point of evaluation, a vector with length 1 |
| cs | The coefficients stored in the order \( [(0, 0), (1, -1), (1, 0), (1, 1) (2, -2), (2, -1), ... (n, n)] \) |
| deg | The degree |
Definition at line 213 of file Sphericals.hpp.
| double Bembel::evaluate_solid_sphericals | ( | Eigen::Vector3d | x, |
| Eigen::VectorXd | cs, | ||
| unsigned int | deg, | ||
| bool | grad | ||
| ) |
Evaluates the series \( \sum_{n = 0}^{\rm deg} \sum_{m = -n}^n |x|^n c_m^n Y_m^n( \frac{x}{|x|}) \) for real coefficients cs if grad is false, and \( \sum_{n = 0}^{\rm deg} \sum_{m = -n}^n n |x|^{n-1} c_m^n Y_m^n(\frac{x}{|x|}) \) for real coefficients cs if grad is true.
| x | The point of evaluation |
| cs | The coefficients stored in the order \( [(0, 0), (1, -1), (1, 0), (1, 1) (2, -2), (2, -1), ... (n, n)] \) |
| deg | The degree |
| grad | distinguish the two cases. |
Definition at line 128 of file Sphericals.hpp.
| double Bembel::evaluate_sphericals | ( | Eigen::Vector3d | x, |
| Eigen::VectorXd | cs, | ||
| unsigned int | deg | ||
| ) |
Evaluates the series \( \sum_{n = 0}^{\rm deg} \sum_{m = -n}^n c_m^n Y_m^n(x) \) for real coefficients, with the convenction that \( Y_{-m}^n := \overline{Y_m^n} \).
| x | The point of evaluation, a vector with length 1 |
| cs | The coefficients stored in the order \( [(0, 0), (1, -1), (1, 0), (1, 1) (2, -2), (2, -1), ... (n, n)] \) |
| deg | The degree |
Definition at line 68 of file Sphericals.hpp.
| Eigen::VectorXd Bembel::getCoefficients | ( | double | precision | ) |
Calculates the coefficients for the solid harmonics expansion of the periodic kernel \( k(x) = \sum_{m \in \{-1, 0, 1\}^3} \frac{1}{4 \pi |x-m|} + \frac{|x|^2}{6} + \sum_{n = 0}^{N} \sum_{m = -n}^n c_m^n \phi_n^m(x) \).
| precision | The desired mean error from periodicity. |
Definition at line 60 of file Coefficients.hpp.
|
inline |
Returns the degree of the sphericals expansion given a precision. Can be extended, use even numbers only!
Definition at line 221 of file Coefficients.hpp.
| Eigen::VectorXd Bembel::getDisplacement | ( | Eigen::MatrixXd | ps_l, |
| Eigen::MatrixXd | ps_f, | ||
| Eigen::MatrixXd | ps_b | ||
| ) |
Returns the right-hand side for the homogenised coefficient calculation.
| ps_l | Points on the left side, i.e. \( \subseteq \{-0.5\} \times [-0.5, 0.5] \times [-0.5, 0.5] \). |
| ps_f | Points on the front side, i.e. \( \subseteq [-0.5, 0.5] \times \{-0.5\} \times [-0.5, 0.5] \). |
| ps_b | Points on the bottom side, i.e. \( \subseteq [-0.5, 0.5] \times [-0.5, 0.5] \times \{-0.5\} \). |
Definition at line 238 of file Coefficients.hpp.
|
constexpr |
Definition at line 75 of file DifferentialFormEnum.hpp.
|
constexpr |
Definition at line 67 of file DifferentialFormEnum.hpp.
|
inlinenoexcept |
method to generate textfile for the geometry
| file_name | name of new textfile |
| number_of_patches | overall number of patches |
Definition at line 130 of file GeometryIO.hpp.
| Eigen::Matrix< double, Eigen::Dynamic, 2 > Bembel::spherical_harmonics_full | ( | Eigen::Vector3d | x, |
| unsigned int | N | ||
| ) |
Calculates the the spherical harmonics \( Y_n^m(\frac{x}{|x|}) \), ordered by \( [Y_0^0, \, Y_1^{-1}, \, Y_1^0, \, Y_1^1, \, Y_2^{-2}, ..., Y_N^N] \).
| x | the point of evaluation |
| N | the maximal degree |
Definition at line 388 of file Sphericals.hpp.
|
noexcept |
method to write Patch information into textfile
| knt1 | knotVector1 |
| knt2 | knotVector2 |
| xyzw | Vector with x,y,z,w Matices |
| file_name | filename |
| current_patch_number | current patch number |
Definition at line 154 of file GeometryIO.hpp.
| std::function<double(const Eigen::Vector2d &, const Eigen::Vector2d &)> Bembel::DummyOperator_test_function |
Definition at line 34 of file DummyOperator.hpp.