46 template<
unsigned int degree,
unsigned int dim>
63 const double basis_value(
unsigned int i,
const arma::vec::fixed<dim> &p)
const;
70 const arma::vec::fixed<dim>
basis_grad(
unsigned int i,
const arma::vec::fixed<dim> &p)
const;
100 template<
unsigned int degree,
unsigned int dim>
155 template <
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
175 double basis_value(
const unsigned int i,
const arma::vec::fixed<dim> &p)
const;
182 arma::vec::fixed<dim>
basis_grad(
const unsigned int i,
const arma::vec::fixed<dim> &p)
const;
187 arma::vec::fixed<dim>
basis_vector(
const unsigned int i,
const arma::vec::fixed<dim> &p)
const;
192 arma::mat::fixed<dim,dim>
basis_grad_vector(
const unsigned int i,
const arma::vec::fixed<dim> &p)
const;
212 template <
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
233 double basis_value(
const unsigned int i,
const arma::vec::fixed<dim> &p)
const;
240 arma::vec::fixed<dim>
basis_grad(
const unsigned int i,
const arma::vec::fixed<dim> &p)
const;
245 arma::vec::fixed<dim>
basis_vector(
const unsigned int i,
const arma::vec::fixed<dim> &p)
const;
250 arma::mat::fixed<dim,dim>
basis_grad_vector(
const unsigned int i,
const arma::vec::fixed<dim> &p)
const;
266 template<
unsigned int degree,
unsigned int dim>
277 arma::uvec::fixed<dim> pows;
280 unsigned int degree_sum=0;
285 powers.push_back(pows);
288 for(i_dim=0; i_dim < dim; i_dim++) {
289 if (degree_sum < degree) {
294 degree_sum-=pows[i_dim];
298 if (i_dim == dim)
break;
302 template<
unsigned int degree,
unsigned int dim>
305 ASSERT(i<=powers.size(),
"Index of basis function is out of range.");
309 for (
unsigned int j=0; j<dim; j++)
310 v *= pow(p[j], (
int) powers[i][j]);
316 template<
unsigned int degree,
unsigned int dim>
319 ASSERT(i<=powers.size(),
"Index of basis function is out of range.");
321 arma::vec::fixed<dim> grad;
323 for (
unsigned int j=0; j<dim; j++)
325 grad[j] = powers[i][j];
326 if (powers[i][j] == 0)
continue;
328 for (
unsigned int k=0; k<dim; k++)
330 grad[j] *= pow(p[k], (
int) (k==j?powers[i][k]-1:powers[i][k]));
346 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
351 for (
int i=0; i<=dim; i++)
353 number_of_dofs += dof_distribution.number_of_single_dofs[i]
354 +2*dof_distribution.number_of_pairs[i]
355 +3*dof_distribution.number_of_triples[i]
356 +6*dof_distribution.number_of_sextuples[i];
358 number_of_single_dofs[i] = dof_distribution.number_of_single_dofs[i];
359 number_of_pairs[i] = dof_distribution.number_of_pairs[i];
360 number_of_triples[i] = dof_distribution.number_of_triples[i];
361 number_of_sextuples[i] = dof_distribution.number_of_sextuples[i];
364 for (
int i=0; i<dof_distribution.unit_support_points.size(); i++)
365 unit_support_points.push_back(dof_distribution.unit_support_points[i]);
369 this->compute_node_matrix();
372 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
375 ASSERT(i <= number_of_dofs,
"Index of basis function is out of range.");
376 return poly_space.basis_value(i, p);
379 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
382 ASSERT(i <= number_of_dofs,
"Index of basis function is out of range.");
383 return poly_space.basis_grad(i, p);
386 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
389 ASSERT(
false,
"basis_vector() may not be called for scalar finite element.");
392 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
395 ASSERT(
false,
"basis_grad_vector() may not be called for scalar finite element.");
398 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
417 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
422 number_of_dofs += dof_distribution.number_of_dofs;
424 number_of_single_dofs[dim] = number_of_dofs;
426 for (
unsigned int i=0; i<dof_distribution.unit_support_points.size(); i++)
427 unit_support_points.push_back(dof_distribution.unit_support_points[i]);
431 this->compute_node_matrix();
434 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
437 ASSERT(i <= number_of_dofs,
"Index of basis function is out of range.");
438 return poly_space.basis_value(i, p);
441 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
444 ASSERT(i <= number_of_dofs,
"Index of basis function is out of range.");
445 return poly_space.basis_grad(i, p);
448 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
451 ASSERT(
false,
"basis_vector() may not be called for scalar finite element.");
452 return arma::vec::fixed<dim>();
455 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
458 ASSERT(
false,
"basis_grad_vector() may not be called for scalar finite element.");
459 return arma::mat::fixed<dim,dim>();
462 template<
unsigned int degree,
unsigned int dim,
unsigned int spacedim>
478 #ifndef DOXYGEN_SHOULD_SKIP_THIS
unsigned int number_of_single_dofs[dim+1]
Number of single dofs at one geometrical entity of the given dimension (point, line, triangle, tetrahedron).
arma::vec::fixed< dim > basis_grad(const unsigned int i, const arma::vec::fixed< dim > &p) const
Returns the gradient of the ith basis function at the point p.
Space of polynomial functions.
unsigned int number_of_dofs
Total number of degrees of freedom at one finite element.
PolynomialSpace()
Constructor.
double basis_value(const unsigned int i, const arma::vec::fixed< dim > &p) const
Returns the ith basis function evaluated at the point p.
unsigned int number_of_triples[dim+1]
Number of triples of dofs associated to one triangle.
unsigned int number_of_sextuples[dim+1]
Number of sextuples of dofs associated to one triangle.
unsigned int number_of_pairs[dim+1]
Number of pairs of dofs at one geometrical entity of the given dimension (applicable to lines and tri...
Discontinuous Lagrangean finite element on dim dimensional simplex.
DofDistribution< degree, dim > dof_distribution
The auxiliary dof distribution.
Global macros to enhance readability and debugging, general constants.
PolynomialSpace< degree, dim > poly_space
The auxiliary polynomial space.
arma::vec::fixed< dim > basis_vector(const unsigned int i, const arma::vec::fixed< dim > &p) const
The vector variant of basis_value must be implemented but may not be used.
virtual ~FE_P_disc()
Destructor.
DofDistribution()
Constructor.
double basis_value(const unsigned int i, const arma::vec::fixed< dim > &p) const
Returns the ith basis function evaluated at the point p.
const arma::vec::fixed< dim > basis_grad(unsigned int i, const arma::vec::fixed< dim > &p) const
Gradient of the i th basis function at point p.
Conforming Lagrangean finite element on dim dimensional simplex.
DofDistribution< degree, dim > dof_distribution
The auxiliary dof distribution.
const double basis_value(unsigned int i, const arma::vec::fixed< dim > &p) const
Value of the i th basis function at point p.
arma::mat::fixed< dim, dim > basis_grad_vector(const unsigned int i, const arma::vec::fixed< dim > &p) const
The vector variant of basis_grad must be implemented but may not be used.
arma::mat::fixed< dim, dim > basis_grad_vector(const unsigned int i, const arma::vec::fixed< dim > &p) const
The vector variant of basis_grad must be implemented but may not be used.
Distribution of dofs for polynomial finite elements.
Abstract class for description of finite elements.
arma::vec::fixed< dim > basis_grad(const unsigned int i, const arma::vec::fixed< dim > &p) const
Returns the gradient of the ith basis function at the point p.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
arma::vec::fixed< dim > basis_vector(const unsigned int i, const arma::vec::fixed< dim > &p) const
The vector variant of basis_value must be implemented but may not be used.
std::vector< arma::uvec::fixed< dim > > powers
Coefficients of basis functions.
std::vector< arma::vec::fixed< dim > > unit_support_points
Support points.
PolynomialSpace< degree, dim > poly_space
The auxiliary polynomial space.