LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
point.cc
1#include "point.h"
2
3namespace lf::geometry {
4Eigen::MatrixXd Point::Global(const Eigen::MatrixXd& local) const {
5 LF_ASSERT_MSG(local.rows() == 0, "local.rows() != 0");
6 return coord_.replicate(1, local.cols());
7}
8
9// NOLINTNEXTLINE(misc-unused-parameters)
10Eigen::MatrixXd Point::Jacobian(const Eigen::MatrixXd& local) const {
11 return Eigen::MatrixXd::Zero(DimGlobal(), 0);
12}
13
15 const Eigen::MatrixXd& local) const { // NOLINT(misc-unused-parameters)
16 LF_VERIFY_MSG(false, "JacobianInverseGramian undefined for points.");
17}
18
19Eigen::VectorXd Point::IntegrationElement(const Eigen::MatrixXd& local) const {
20 return Eigen::VectorXd::Ones(local.cols());
21}
22
23std::unique_ptr<Geometry> Point::SubGeometry(dim_t codim, dim_t i) const {
24 if (codim == 0 && i == 0) {
25 return std::make_unique<Point>(coord_);
26 }
27 LF_VERIFY_MSG(false, "codim or i out of bounds.");
28}
29
30std::vector<std::unique_ptr<Geometry>> Point::ChildGeometry(
31 const RefinementPattern& ref_pattern, lf::base::dim_t codim) const {
32 LF_VERIFY_MSG(codim == 0, "Only codim = 0 allowed for a point");
33 std::vector<std::unique_ptr<Geometry>> child_geo_uptrs{};
34 LF_VERIFY_MSG(ref_pattern.RefEl() == lf::base::RefEl::kPoint(),
35 "ref_patern.RefEl() = " << ref_pattern.RefEl().ToString());
36 LF_VERIFY_MSG(ref_pattern.NumChildren(0) == 1,
37 "ref_pattern.NumChildren() = " << ref_pattern.NumChildren(0));
38 // The only way to "refine" a point is to copy it
39 child_geo_uptrs.push_back(std::make_unique<Point>(coord_));
40 return child_geo_uptrs;
41}
42
43} // namespace lf::geometry
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
std::string ToString() const
Return a string representation of this Reference element.
Definition ref_el.h:458
base::RefEl::dim_t dim_t
std::unique_ptr< Geometry > SubGeometry(dim_t codim, dim_t i) const override
Construct a new Geometry() object that describes the geometry of the i-th sub-entity with codimension...
Definition point.cc:23
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition point.cc:10
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition point.cc:4
Eigen::VectorXd coord_
Definition point.h:41
std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pattern, base::dim_t codim) const override
the child geometry is just a copy of the point geometry
Definition point.cc:30
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition point.h:16
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition point.cc:14
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Definition point.cc:19
Abstract interface class for encoding topological local refinement
virtual lf::base::size_type NumChildren(lf::base::dim_t codim) const =0
provide number of child entities of a given co-dimension to be created by refinement
lf::base::RefEl RefEl() const
Returns topological type of entity for which the current object is set up.
unsigned int dim_t
type for dimensions and co-dimensions and numbers derived from them
Definition types.h:32
Defines the Geometry interface and provides a number of classes that implement this interface + addit...
Definition compose.h:13