45 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesData;
57 template<>
inline double determinant(
const arma::mat::fixed<1,2> &M)
59 return sqrt(M(0,0)*M(0,0)+M(0,1)*M(0,1));
62 template<>
inline double determinant(
const arma::mat::fixed<2,1> &M)
64 return sqrt(M(0,0)*M(0,0)+M(1,0)*M(1,0));
67 template<>
inline double determinant(
const arma::mat::fixed<0,3> &M)
72 template<>
inline double determinant(
const arma::mat::fixed<3,0> &M)
77 template<>
inline double determinant(
const arma::mat::fixed<1,3> &M)
79 return sqrt(M(0,0)*M(0,0)+M(0,1)*M(0,1)+M(0,2)*M(0,2));
82 template<>
inline double determinant(
const arma::mat::fixed<3,1> &M)
84 return sqrt(M(0,0)*M(0,0)+M(1,0)*M(1,0)+M(2,0)*M(2,0));
87 template<>
inline double determinant(
const arma::mat::fixed<2,3> &M)
89 return sqrt((M(0,0)*M(0,0)+M(0,1)*M(0,1)+M(0,2)*M(0,2))*(M(1,0)*M(1,0)+M(1,1)*M(1,1)+M(1,2)*M(1,2))
90 -(M(0,0)*M(1,0)+M(0,1)*M(1,1)+M(0,2)*M(1,2))*(M(0,0)*M(1,0)+M(0,1)*M(1,1)+M(0,2)*M(1,2)));
93 template<>
inline double determinant(
const arma::mat::fixed<3,2> &M)
95 return sqrt((M(0,0)*M(0,0)+M(1,0)*M(1,0)+M(2,0)*M(2,0))*(M(0,1)*M(0,1)+M(1,1)*M(1,1)+M(2,1)*M(2,1))
96 -(M(0,0)*M(0,1)+M(1,0)*M(1,1)+M(2,0)*M(2,1))*(M(0,0)*M(0,1)+M(1,0)*M(1,1)+M(2,0)*M(2,1)));
99 template<
unsigned int n>
inline double determinant(
const arma::mat::fixed<n,n> &M)
128 template<
unsigned int dim,
unsigned int spacedim>
200 template<
unsigned int dim,
unsigned int spacedim>
inline
211 arma::vec::fixed<dim+1> el_bar_coords;
212 arma::vec::fixed<dim> side_bar_coords;
214 for (
unsigned int k=0; k<q.
size(); k++)
218 el_bar_coords.zeros();
221 for (
int j=0; j<dim-1; j++)
223 side_bar_coords(j) = (subq.
point(k))[j];
224 lambda += (subq.
point(k))[j];
226 side_bar_coords(dim-1) = 1.-lambda;
229 for (
unsigned int i=0; i<dim; i++)
231 q.
set_point(k, el_bar_coords.subvec(0,dim-1));
void set_weight(const unsigned int i, const double w)
Sets individual quadrature weight.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
virtual void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &cell, unsigned int sid, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)=0
Calculates the mapping data related to a given side, namely the jacobian determinants and the normal ...
virtual ~Mapping()
Destructor.
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
virtual MappingInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)=0
Calculates the mapping data on the reference cell.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Base class for quadrature rules on simplices in arbitrary dimensions.
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
double weight(const unsigned int i) const
Returns the ith weight.
Abstract class for the mapping between reference and actual cell.
void resize(const unsigned int n_q_points)
Modify the number of quadrature points.
void transform_subquadrature(unsigned int sid, unsigned int pid, const Quadrature< dim-1 > &subq, Quadrature< dim > &q)
Creates a cell dim-dimensional quadrature from side (dim-1)-dimensional quadrature.
std::vector< arma::vec > bar_coords
Auxiliary array of barycentric coordinates of quadrature points.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
void set_point(const unsigned int i, const arma::vec::fixed< dim > &p)
Sets individual quadrature point coordinates.
const unsigned int size() const
Returns number of quadrature points.
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
virtual void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)=0
Calculates the mapping data and stores them in the provided structures.
virtual UpdateFlags update_each(UpdateFlags flags)=0
Decides which additional quantities have to be computed for each cell.
Mapping data that can be precomputed on the actual cell.