Bembel
Bembel Namespace Reference

Routines for the evalutation of pointwise errors. More...

Detailed Description

Routines for the evalutation of pointwise errors.

Function Documentation

◆ calculateFirstCoefficient()

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.

Parameters
csThe already calculated other coefficients
degThe degree
ps_lGauss quadrature points on the left side of the cube
ps_fGauss quadrature points on the front side of the cube
ps_bGauss quadrature points on the bottom side of the cube

Definition at line 346 of file Coefficients.hpp.

◆ double_to_string()

std::string Bembel::double_to_string ( double  d,
const int  precision 
)

Convert a double to a string with given precision.

Parameters
ddouble to be converted.
precisionPrecision of the conversion.

Definition at line 253 of file GeometryIGS.hpp.

◆ evaluate_dsolid_sphericals()

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.

Parameters
xThe point of evaluation
csThe coefficients stored in the order \( [(0, 0), (1, -1), (1, 0), (1, 1) (2, -2), (2, -1), ... (n, n)] \)
degThe degree

Definition at line 302 of file Sphericals.hpp.

◆ evaluate_dsphericals()

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.

Parameters
xThe point of evaluation, a vector with length 1
csThe coefficients stored in the order \( [(0, 0), (1, -1), (1, 0), (1, 1) (2, -2), (2, -1), ... (n, n)] \)
degThe degree

Definition at line 213 of file Sphericals.hpp.

◆ evaluate_solid_sphericals()

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.

Parameters
xThe point of evaluation
csThe coefficients stored in the order \( [(0, 0), (1, -1), (1, 0), (1, 1) (2, -2), (2, -1), ... (n, n)] \)
degThe degree
graddistinguish the two cases.

Definition at line 128 of file Sphericals.hpp.

◆ evaluate_sphericals()

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} \).

See also
K. Giebermann. Schnelle Summationsverfahren zur numerischen Lösung von Integralgleichungen für Streuprobleme im \( \mathbb{R}^3 \). PhD thesis, Universität Karlsruhe (TH), Germany, 1997.
R. von Rickenbach. Boundary Element Methods for Shape Optimisation in Homogenisation. MSc thesis, Universität Basel, Switzerland, 2021.
Parameters
xThe point of evaluation, a vector with length 1
csThe coefficients stored in the order \( [(0, 0), (1, -1), (1, 0), (1, 1) (2, -2), (2, -1), ... (n, n)] \)
degThe degree

Definition at line 68 of file Sphericals.hpp.

◆ getCoefficients()

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) \).

Parameters
precisionThe desired mean error from periodicity.
See also
A. Barnett and L. Greengard. A new integral representation for quasi-periodic fields and its application to two-dimensional band structure calculations. Journal of Computational Physics, 229(19):6898–6914, 2010.
P. Cazeaux and O. Zahm. A fast boundary element method for the solution of periodic many-inclusion problems via hierarchical matrix techniques. ESAIM: Proceedings and Surveys, 48:156–168, 2015.
R. von Rickenbach. Boundary Element Methods for Shape Optimisation in Homogenisation. MSc thesis, Universität Basel, 2021.

Definition at line 60 of file Coefficients.hpp.

◆ getDegree()

unsigned int Bembel::getDegree ( double  precision)
inline

Returns the degree of the sphericals expansion given a precision. Can be extended, use even numbers only!

See also
P. Cazeaux and O. Zahm. A fast boundary element method for the solution of periodic many-inclusion problems via hierarchical matrix techniques. ESAIM: Proceedings and Surveys, 48:156–168, 2015.
R. von Rickenbach. Boundary Element Methods for Shape Optimisation in Homogenisation. MSc thesis, Universität Basel, 2021.

Definition at line 221 of file Coefficients.hpp.

◆ getDisplacement()

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.

Parameters
ps_lPoints on the left side, i.e. \( \subseteq \{-0.5\} \times [-0.5, 0.5] \times [-0.5, 0.5] \).
ps_fPoints on the front side, i.e. \( \subseteq [-0.5, 0.5] \times \{-0.5\} \times [-0.5, 0.5] \).
ps_bPoints 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.

◆ getFunctionSpaceOutputDimension()

template<unsigned int DF>
constexpr int Bembel::getFunctionSpaceOutputDimension ( )
constexpr
Deprecated:
Use DifferentialFormTraits instead.

Definition at line 75 of file DifferentialFormEnum.hpp.

◆ getFunctionSpaceVectorDimension()

template<unsigned int DF>
constexpr int Bembel::getFunctionSpaceVectorDimension ( )
constexpr
Deprecated:
Use DifferentialFormTraits instead.

Definition at line 67 of file DifferentialFormEnum.hpp.

◆ MakeFile()

void Bembel::MakeFile ( const std::string &  file_name,
int  number_of_patches 
)
inlinenoexcept

method to generate textfile for the geometry

Parameters
file_namename of new textfile
number_of_patchesoverall number of patches

Definition at line 130 of file GeometryIO.hpp.

◆ PatchShredder() [1/2]

std::vector<Patch> Bembel::PatchShredder ( const Patch patch)
inlinenoexcept

This function cuts a patch along internal knots, if any.

Parameters
patchPatch which gets checked and cut if there are internal knots.
Returns
Vector of patches if it got cut, otherwise the vector is of size 1.

Definition at line 386 of file Patch.hpp.

◆ PatchShredder() [2/2]

std::vector<Patch> Bembel::PatchShredder ( const std::vector< Patch > &  patches)
inlinenoexcept

This function cuts all patches along internal knots, if any.

Parameters
patchesVector with patches.
Returns
Vector of patches if it got cut.

Definition at line 437 of file Patch.hpp.

◆ spherical_harmonics_full()

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] \).

Parameters
xthe point of evaluation
Nthe maximal degree

Definition at line 388 of file Sphericals.hpp.

◆ WritePatch()

void Bembel::WritePatch ( const std::string &  file_name,
int  current_patch_number,
const std::vector< Eigen::MatrixXd > &  xyzw,
const std::vector< double > &  knt1,
const std::vector< double > &  knt2 
)
noexcept

method to write Patch information into textfile

Parameters
knt1knotVector1
knt2knotVector2
xyzwVector with x,y,z,w Matices
file_namefilename
current_patch_numbercurrent patch number

Definition at line 154 of file GeometryIO.hpp.

Variable Documentation

◆ DummyOperator_test_function

std::function<double(const Eigen::Vector2d &, const Eigen::Vector2d &)> Bembel::DummyOperator_test_function
Initial value:
=
[](const Eigen::Vector2d &x, const Eigen::Vector2d &y) { return 1.; }

Definition at line 34 of file DummyOperator.hpp.