44 template<
unsigned int dim,
unsigned int spacedim>
50 template<
unsigned int dim,
unsigned int spacedim>
55 for (
unsigned int i = 0; i <= dim; i++)
57 number_of_single_dofs[i] = 0;
58 number_of_pairs[i] = 0;
59 number_of_triples[i] = 0;
60 number_of_sextuples[i] = 0;
64 template<
unsigned int dim,
unsigned int spacedim>
inline
67 return number_of_dofs;
70 template<
unsigned int dim,
unsigned int spacedim>
inline
74 ASSERT(object_dim >= 0 && object_dim <= dim,
75 "Object type number is out of range.");
79 return number_of_single_dofs[object_dim];
81 return number_of_pairs[object_dim];
83 return number_of_triples[object_dim];
85 return number_of_sextuples[object_dim];
91 template<
unsigned int dim,
unsigned int spacedim>
inline
94 ASSERT_EQUAL(get_generalized_support_points().size(), number_of_dofs);
96 arma::mat M(number_of_dofs, number_of_dofs);
98 for (
unsigned int i = 0; i < number_of_dofs; i++)
99 for (
unsigned int j = 0; j < number_of_dofs; j++) {
100 M(j, i) = basis_value(j, get_generalized_support_points()[i]);
103 node_matrix = arma::inv(M);
106 template<
unsigned int dim,
unsigned int spacedim>
113 arma::vec values(number_of_dofs);
115 for (
unsigned int i=0; i<q.
size(); i++)
117 for (
unsigned int j=0; j<number_of_dofs; j++)
118 values[j] = basis_value(j, q.
point(i));
125 arma::mat grads(number_of_dofs, dim);
127 for (
unsigned int i=0; i<q.
size(); i++)
129 for (
unsigned int j=0; j<number_of_dofs; j++)
130 grads.row(j) = arma::trans(basis_grad(j, q.
point(i)));
138 template<
unsigned int dim,
unsigned int spacedim>
inline
149 template<
unsigned int dim,
unsigned int spacedim>
inline
158 for (
unsigned int i = 0; i < q.
size(); i++)
165 for (
unsigned int i = 0; i < q.
size(); i++)
172 template<
unsigned int dim,
unsigned int spacedim>
175 if (generalized_support_points.size() > 0)
177 return generalized_support_points;
181 return unit_support_points;
186 template<
unsigned int dim,
unsigned int spacedim>
std::vector< arma::vec > shape_values
Shape functions evaluated at the quadrature points.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
const std::vector< arma::vec::fixed< dim > > & get_generalized_support_points()
Returns either the generalized support points (if they are defined) or the unit support points...
void init()
Clears all internal structures.
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Class FEValues calculates finite element data on the actual cells such as shape function values...
std::vector< arma::mat > basis_grads
Precomputed gradients of basis functions at the quadrature points.
Base class for quadrature rules on simplices in arbitrary dimensions.
virtual UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
virtual ~FiniteElement()
Destructor.
#define ASSERT_EQUAL(a, b)
virtual void compute_node_matrix()
Initializes the node_matrix for computing the coefficients of the raw basis functions from values at ...
std::vector< arma::mat > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
Basic definitions of numerical quadrature rules.
Shape function gradients.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
const unsigned int n_object_dofs(unsigned int object_dim, DofMultiplicity multiplicity)
Returns the number of single dofs/dof pairs/triples/sextuples that lie on a single geometric entity o...
const unsigned int size() const
Returns number of quadrature points.
std::vector< arma::vec > basis_values
Precomputed values of basis functions at the quadrature points.
DofMultiplicity
Multiplicity of finite element dofs.
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
FiniteElement()
Constructor.
const unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
virtual FEInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)
Calculates the data on the reference cell.
Structure for storing the precomputed finite element data.
Abstract class for description of finite elements.
virtual void fill_fe_values(const Quadrature< dim > &q, FEInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure...
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...