44 template<
unsigned int dim,
unsigned int spacedim>
49 template<
unsigned int dim,
unsigned int spacedim>
65 for (
unsigned int i=0; i<dim; i++)
75 vec::fixed<dim+1> basis;
77 for (
unsigned int i=0; i<q.
size(); i++)
80 for (
unsigned int j=0; j<dim; j++)
82 basis[0] -= q.
point(i)[j];
83 basis[j+1] = q.
point(i)[j];
94 template<
unsigned int dim,
unsigned int spacedim>
111 template<
unsigned int dim,
unsigned int spacedim>
117 mat::fixed<spacedim,dim+1> coords;
118 mat::fixed<spacedim,dim> jac;
127 for (
unsigned int n=0; n<dim+1; n++)
128 for (
unsigned int c=0; c<spacedim; c++)
129 coords(c,n) = cell->node[n]->point()[c];
142 for (
unsigned int i=0; i<q.
size(); i++)
152 for (
unsigned int i=0; i<q.
size(); i++)
157 for (
unsigned int i=0; i<q.
size(); i++)
164 mat::fixed<dim,spacedim> ijac;
173 for (
unsigned int i=0; i<q.
size(); i++)
181 vec::fixed<dim+1> basis;
182 for (
unsigned int i=0; i<q.
size(); i++)
187 template<
unsigned int dim,
unsigned int spacedim>
194 mat::fixed<spacedim,dim+1> coords;
203 for (
unsigned int n=0; n<dim+1; n++)
204 for (
unsigned int c=0; c<spacedim; c++)
205 coords(c,n) = cell->node[n]->point()[c];
214 mat::fixed<spacedim,dim> jac = coords*grad;
218 for (
unsigned int i=0; i<q.
size(); i++)
225 for (
unsigned int i=0; i<q.
size(); i++)
232 mat::fixed<dim,spacedim> ijac;
241 for (
unsigned int i=0; i<q.
size(); i++)
247 vec::fixed<spacedim> n_cell;
249 n_cell = n_cell/norm(n_cell,2);
250 for (
unsigned int i=0; i<q.
size(); i++)
260 for (
unsigned int i=0; i<q.
size(); i++)
273 mat::fixed<spacedim,dim> side_coords;
274 mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac;
278 for (
unsigned int n=0; n<dim; n++)
279 for (
unsigned int c=0; c<spacedim; c++)
280 side_coords(c,n) = cell->side(sid)->node(n)->point()[c];
281 side_jac = side_coords * grad.submat(0,0,dim-1,dim-2);
286 for (
unsigned int i=0; i<q.
size(); i++)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Transformed quadrature weight for cell sides.
MappingInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)
Initializes the structures and computes static data.
Determinant of the Jacobian.
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.
UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Class FEValues calculates finite element data on the actual cells such as shape function values...
std::vector< double > JxW_values
Transformed quadrature weights.
Base class for quadrature rules on simplices in arbitrary dimensions.
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Normal vectors to the element at the quadrature points lying on a side.
void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Calculates the mapping data on the actual cell.
Transformed quadrature points.
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.
double weight(const unsigned int i) const
Returns the ith weight.
void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &cell, unsigned int sid, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Calculates the mapping data on a side of a cell.
std::vector< arma::vec > bar_coords
Auxiliary array of barycentric coordinates of quadrature points.
Basic definitions of numerical quadrature rules.
Affine mapping between reference and actual cell.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
const unsigned int size() const
Returns number of quadrature points.
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Mapping data that can be precomputed on the actual cell.
static arma::vec::fixed< dim > normal_vector(unsigned int sid)
Transformed quadrature weights.