Flow123d
release_3.0.0-973-g92f55e826
|
Go to the documentation of this file.
40 unsigned int degree_sum=0;
47 for(i_dim=0; i_dim <
dim; i_dim++) {
48 if (degree_sum < degree) {
53 degree_sum-=pows[i_dim];
57 if (i_dim ==
dim)
break;
63 const arma::vec &point,
64 unsigned int comp_index
72 for (
unsigned int j=0; j<this->
space_dim_; j++)
73 v *= pow(point[j], (
int)
powers[i][j]);
81 unsigned int comp_index
89 for (
unsigned int j=0; j<this->
space_dim_; j++)
92 if (
powers[i][j] == 0)
continue;
94 for (
unsigned int k=0; k<this->
space_dim_; k++)
96 grad[j] *= pow(p[k], (
int) (k==j?
powers[i][k]-1:
powers[i][k]));
112 template<
unsigned int dim>
115 if (degree_ == 0 ||
dim == 0)
120 arma::vec coords = arma::ones<arma::vec>(
dim+1)/(
dim+1);
121 this->dofs_.push_back(
Dof(
dim, 0, coords, { 1 },
Value));
131 arma::uvec ubc = arma::zeros<arma::uvec>(
dim+1);
148 while (ubc[c] == 0) c++;
150 if (c ==
dim) finish =
true;
160 for (
auto ubc : uvbc)
167 arma::uvec nonzeros = find(ubc);
169 arma::vec coords = arma::conv_to<arma::vec>::from(ubc);
174 unsigned int n_face_idx;
175 switch (
dim-zeros.first) {
189 this->dofs_.push_back(
Dof(nonzeros.size()-1, n_face_idx, coords, { 1 },
Value));
197 template<
unsigned int dim>
225 template<
unsigned int dim>
230 for (
unsigned int i=0; i<this->
dofs_.size(); i++)
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
FE_P_disc(unsigned int degree)
Constructor.
static std::pair< unsigned int, unsigned int > zeros_positions(const BaryPoint &barycentric, double tolerance=std::numeric_limits< double >::epsilon() *2)
std::vector< Dof > dofs_
Set of degrees of freedom (functionals) defining the FE.
unsigned int dim() const
Return dimension of element appropriate to the side.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
virtual void compute_node_matrix()
Initializes the node_matrix for computing the coefficients of the shape functions in the raw basis of...
Conforming Lagrangean finite element on dim dimensional simplex.
const unsigned int dim() const override
Dimension of function space (number of basis functions).
FE_P(unsigned int degree)
Constructor.
const arma::vec basis_grad(unsigned int basis_index, const arma::vec &point, unsigned int comp_index=0) const override
Gradient of the i th basis function at point point.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
std::vector< arma::uvec > powers
Coefficients of basis functions.
PolynomialSpace(unsigned int degree, unsigned int dim)
Constructor.
std::shared_ptr< FunctionSpace > function_space_
Function space defining the FE.
unsigned int space_dim_
Space dimension of function arguments (i.e. 1, 2 or 3).
const double basis_value(unsigned int basis_index, const arma::vec &point, unsigned int comp_index=0) const override
Value of the i th basis function at point point.
void setup_components()
Initialize vectors with information about components of basis functions.