LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
ref_el.cc
1#include "ref_el.h"
2
3namespace lf::base {
4
5const Eigen::MatrixXd RefEl::ncoords_point_dynamic_ = Eigen::VectorXd(0);
6
7const Eigen::MatrixXd RefEl::ncoords_segment_dynamic_ =
8 Eigen::Vector2d(0, 1).transpose();
9
10const Eigen::MatrixXd RefEl::ncoords_tria_dynamic_ =
11 (Eigen::MatrixXd{2, 3} << 0, 1, 0, 0, 0, 1).finished();
12
13const Eigen::MatrixXd RefEl::ncoords_quad_dynamic_ =
14 (Eigen::MatrixXd{2, 4} << 0, 1, 1, 0, 0, 0, 1, 1).finished();
15
16// Print function
17void PrintInfo(std::ostream &stream, const RefEl &ref_el, int output_ctrl) {
18 const base::dim_t dim_ref_el = ref_el.Dimension();
19 const base::dim_t no_nodes = ref_el.NumNodes();
20 stream << "Type of reference element: " << ref_el.ToString() << '\n';
21 stream << "Dimension: " << dim_ref_el << '\n';
22 stream << "Number of nodes: " << no_nodes << '\n';
23
24 if (output_ctrl > 0) {
25 // Loop over dimensions
26 for (base::dim_t co_dim = dim_ref_el; co_dim > 0; co_dim--) {
27 base::dim_t num_sub_ent = ref_el.NumSubEntities(co_dim);
28 stream << "Codimension " << co_dim << " has " << num_sub_ent
29 << " entities of type " << ref_el.SubType(co_dim, 0).ToString()
30 << '\n';
31
32 if (output_ctrl > 10) {
33 for (; num_sub_ent > 0; num_sub_ent--) {
34 const std::int32_t sub_ent =
35 static_cast<std::int32_t>(num_sub_ent) - 1;
36 stream << " Subentity " << sub_ent << " is of type "
37 << ref_el.SubType(co_dim, 0).ToString();
38
39 if (ref_el.SubType(co_dim, 0) == RefEl::kPoint() &&
40 output_ctrl > 20) {
41 stream << " and has coordinates ["
42 << ref_el.NodeCoords().col(sub_ent)[0] << " "
43 << ref_el.NodeCoords().col(sub_ent)[1] << "]" << '\n';
44 }
45 stream << '\n';
46 }
47 }
48 }
49 }
50} // end void PrintInfo
51
52} // namespace lf::base
Represents a reference element with all its properties.
Definition ref_el.h:109
static const Eigen::MatrixXd ncoords_quad_dynamic_
Definition ref_el.h:114
static const Eigen::MatrixXd ncoords_segment_dynamic_
Definition ref_el.h:112
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
Definition ref_el.h:241
constexpr size_type NumSubEntities(dim_t sub_codim) const
Get the number of sub-entities of this RefEl with the given codimension.
Definition ref_el.h:308
constexpr dim_t Dimension() const
Return the dimension of this reference element.
Definition ref_el.h:204
static const Eigen::MatrixXd ncoords_tria_dynamic_
Definition ref_el.h:113
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
static const Eigen::MatrixXd ncoords_point_dynamic_
Definition ref_el.h:111
constexpr size_type NumNodes() const
The number of nodes of this reference element.
Definition ref_el.h:213
std::string ToString() const
Return a string representation of this Reference element.
Definition ref_el.h:458
constexpr RefEl SubType(dim_t sub_codim, dim_t sub_index) const
Return the RefEl of the sub-entity with codim sub_codim and index sub_index.
Definition ref_el.h:351
unsigned int dim_t
type for dimensions and co-dimensions and numbers derived from them
Definition types.h:32
Contains basic functionality that is used by other parts of LehrFEM++.
Definition base.h:15
void PrintInfo(std::ostream &stream, const RefEl &ref_el, int output_ctrl)
Definition ref_el.cc:17