32 template<
unsigned int dim,
unsigned int spacedim>
37 template<
unsigned int dim,
unsigned int spacedim>
53 for (
unsigned int i=0; i<dim; i++)
63 vec::fixed<dim+1> basis;
65 for (
unsigned int i=0; i<q.
size(); i++)
68 for (
unsigned int j=0; j<dim; j++)
70 basis[0] -= q.
point(i)[j];
71 basis[j+1] = q.
point(i)[j];
82 template<
unsigned int dim,
unsigned int spacedim>
99 template<
unsigned int dim,
unsigned int spacedim>
105 mat::fixed<spacedim,dim+1> coords;
106 mat::fixed<spacedim,dim> jac;
115 for (
unsigned int n=0; n<dim+1; n++)
116 for (
unsigned int c=0; c<spacedim; c++)
117 coords(c,n) = cell->node[n]->point()[c];
130 for (
unsigned int i=0; i<q.
size(); i++)
140 for (
unsigned int i=0; i<q.
size(); i++)
145 for (
unsigned int i=0; i<q.
size(); i++)
152 mat::fixed<dim,spacedim> ijac;
161 for (
unsigned int i=0; i<q.
size(); i++)
169 vec::fixed<dim+1> basis;
170 for (
unsigned int i=0; i<q.
size(); i++)
175 template<
unsigned int dim,
unsigned int spacedim>
182 mat::fixed<spacedim,dim+1> coords;
191 for (
unsigned int n=0; n<dim+1; n++)
192 for (
unsigned int c=0; c<spacedim; c++)
193 coords(c,n) = cell->node[n]->point()[c];
202 mat::fixed<spacedim,dim> jac = coords*grad;
206 for (
unsigned int i=0; i<q.
size(); i++)
213 for (
unsigned int i=0; i<q.
size(); i++)
220 mat::fixed<dim,spacedim> ijac;
229 for (
unsigned int i=0; i<q.
size(); i++)
235 vec::fixed<spacedim> n_cell;
237 n_cell = n_cell/norm(n_cell,2);
238 for (
unsigned int i=0; i<q.
size(); i++)
248 for (
unsigned int i=0; i<q.
size(); i++)
261 mat::fixed<spacedim,dim> side_coords;
262 mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac;
266 for (
unsigned int n=0; n<dim; n++)
267 for (
unsigned int c=0; c<spacedim; c++)
268 side_coords(c,n) = cell->side(sid)->node(n)->point()[c];
269 side_jac = side_coords * grad.submat(0,0,dim-1,dim-2);
274 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...
static LocalPoint normal_vector(unsigned int sid)
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.
Transformed quadrature weights.