Flow123d
JS_before_hm-1716-g9144da4bf
|
Go to the documentation of this file.
71 #ifndef REF_ELEMENT_HH_
72 #define REF_ELEMENT_HH_
138 IdxVector(std::array<unsigned int,Size> data_in);
140 IdxVector(std::initializer_list<unsigned int> data_in);
142 unsigned int operator[](unsigned int idx) const;
146 template<std::
size_t Size>
156 template <
unsigned int OutDim,
unsigned int InDim>
162 template<
unsigned int dim>
234 static const unsigned int n_lines = (
unsigned int)((dim * (dim + 1)) / 2);
317 double tolerance = std::numeric_limits<double>::epsilon()*2);
334 template<
unsigned int subdim>
static unsigned int count();
343 template<
unsigned int subdim>
static arma::mat::fixed<dim+1,subdim+1>
bary_coords(
unsigned int sid);
349 template<
unsigned int subdim>
static BaryPoint interpolate(arma::vec::fixed<subdim+1> coord,
int sub_simplex_idx);
357 double first_theta,
double second_theta,
double theta);
390 template <
template <
unsigned int OutDim,
unsigned int InDim>
class TInteraction,
unsigned int OutDim,
unsigned int InDim>
391 static const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) >
interact( TInteraction<OutDim,InDim> interaction );
396 template<
unsigned int OutDim,
unsigned int InDim>
397 static const IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) >
interact_(
unsigned int index);
473 template<
unsigned int dim>
474 template<
unsigned int subdim>
inline
477 arma::mat::fixed<dim+1,subdim+1> bary_c;
479 for(
unsigned int i = 0; i < subdim+1; i++){
480 unsigned int nid = interact_<0,subdim>(sid)[i];
481 bary_c.col(i).zeros();
482 bary_c.col(i)(nid) = 1;
489 template<
unsigned int dim>
inline
492 ASSERT_LT_DBG(nid, n_nodes).error(
"Node number is out of range!");
494 arma::vec::fixed<dim> p;
504 template<
unsigned int dim>
505 template<
unsigned int subdim>
562 template<
unsigned int dim>
563 template<
unsigned int subdim>
566 for(
unsigned int i=0; i < RefElement<dim>::count<subdim>(); i++){
567 if(zeros_positions == topology_zeros_[subdim][i])
return i;
569 ASSERT(0).error(
"Undefined zero pattern.");
582 return line_nodes_[i];}
587 return line_nodes_[i];}
592 return line_nodes_[i];}
597 return node_lines_[i];}
602 return node_lines_[i];}
607 return node_lines_[i];}
612 return side_nodes_[i];}
617 return node_sides_[i];}
622 return line_sides_[i];}
627 return side_lines_[i];}
629 template<
unsigned int dim>
template<
unsigned int OutDim,
unsigned int InDim>
632 ASSERT(
false)(dim)(OutDim)(InDim)(i).error(
"Not implemented.");
635 return IdxVector< (InDim>OutDim ? InDim+1 : dim-InDim) >();
639 template<
unsigned int dim>
640 template <
template <
unsigned int OutDim,
unsigned int InDim>
class TInteraction,
unsigned int OutDim,
unsigned int InDim>
643 return interact_<OutDim,InDim>(interaction.i_);
static LocalPoint normal_vector(unsigned int sid)
static unsigned int permutation_index(unsigned int p[n_nodes_per_side])
Armor::ArmaVec< double, dim > FaceBaryPoint
static BaryPoint line_barycentric_interpolation(BaryPoint first_coords, BaryPoint second_coords, double first_theta, double second_theta, double theta)
static unsigned int normal_orientation(unsigned int sid)
Armor::ArmaVec< double, dim+1 > BaryPoint
arma::vec::fixed< dim > LocalPoint
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
static const unsigned int n_sides_per_line
Number of sides with one common line. dim == 3.
static arma::mat::fixed< dim+1, subdim+1 > bary_coords(unsigned int sid)
static const std::vector< IdxVector< n_nodes_per_side > > side_nodes_
[n_sides] For given side, returns nodes indices. For dim == 3.
static const unsigned int n_nodes
Number of nodes.
static BaryPoint interpolate(arma::vec::fixed< subdim+1 > coord, int sub_simplex_idx)
static std::pair< unsigned int, unsigned int > zeros_positions(const BaryPoint &barycentric, double tolerance=std::numeric_limits< double >::epsilon() *2)
static unsigned int count()
static const std::vector< IdxVector< n_lines_per_side > > side_lines_
[n_sides] For given side, returns lines indices. For dim == 3.
static const unsigned int n_sides
Number of sides.
static BarycentricUnitVec make_bary_unit_vec()
static CentersList centers_of_subelements(unsigned int sub_dim)
static const std::vector< std::vector< std::vector< unsigned int > > > nodes_of_subelements
static const IdxVector<(InDim >OutDim ? InDim+1 :dim-InDim) > interact_(unsigned int index)
Internal part of the interact function.
static const std::vector< IdxVector< n_sides_per_node > > node_sides_
[n_nodes] For given node, returns sides indices. For dim == 3.
static LocalPoint bary_to_local(const BaryPoint &bp)
Converts from barycentric to local coordinates.
static const unsigned int n_nodes_per_line
Number of nodes in one line.
static const unsigned int n_sides_per_node
Number of sides with one common line.
static BaryPoint clip(const BaryPoint &barycentric)
std::array< unsigned int, Size > IdxVector
static const IdxVector<(n_lines > n_nodes) ? n_lines :n_nodes > topology_zeros_[dim+1]
static unsigned int oposite_node(unsigned int sid)
static const IdxVector<(InDim >OutDim ? InDim+1 :dim-InDim) > interact(TInteraction< OutDim, InDim > interaction)
const typedef std::vector< LocalPoint > & CentersList
DECLARE_EXCEPTION(ExcInvalidPermutation,<< "Side permutation not found.\n")
static FaceBaryPoint barycentric_on_face(const BaryPoint &barycentric, unsigned int i_face)
static const std::vector< IdxVector< n_nodes_per_line > > line_nodes_
[n_lines] For given line, returns its nodes indices.
static const unsigned int n_lines
Number of lines, i.e. object of dimension dim-2 on the boundary of the reference element.
static LocalPoint node_coords(unsigned int nid)
static const unsigned int n_lines_per_side
Number of lines on boundary of one side.
static unsigned int line_between_faces(unsigned int f1, unsigned int f2)
Interaction(unsigned int i)
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
std::vector< BaryPoint > BarycentricUnitVec
static const std::vector< IdxVector< n_sides_per_line > > line_sides_
[n_lines] For given line, returns sides indices. For dim == 3.
static double side_measure(unsigned int sid)
static const unsigned int n_lines_per_node
Number of lines with one common node.
static constexpr unsigned int n_side_permutations
static const unsigned int n_nodes_per_side
Number of nodes on one side.
typename arma::Col< Type >::template fixed< nr > ArmaVec
static const std::vector< std::vector< unsigned int > > side_permutations
static unsigned int topology_idx(unsigned int zeros_positions)
static const std::vector< IdxVector< n_lines_per_node > > node_lines_
[n_nodes] For given node, returns lines indices.