38 template<
unsigned int spacedim>
43 n_points_(_quadrature.size()),
44 n_sides_(_quadrature.dim() == dim ? 0 : dim+1),
46 side_ref_data(n_sides_)
60 template<
unsigned int spacedim>
76 template<
unsigned int spacedim>
100 for (
unsigned int c=0; c<spacedim; c++)
printf(
"%g ", n[c]);
111 template<
unsigned int spacedim>
118 for (
unsigned int i=0; i<q.
size(); i++)
147 template<
unsigned int spacedim>
173 template<
unsigned int spacedim>
179 data(this->n_points_, update_each(_flags), dim)
184 template<
unsigned int spacedim>
188 if ( _quadrature.
dim() == dim )
191 this->ref_data = this->init_ref_data(_quadrature);
193 else if ( _quadrature.
dim() + 1 == dim )
196 for (
unsigned int sid = 0; sid < this->n_sides_; sid++)
215 this->side_ref_data[sid] = this->init_ref_data(side_quad);
221 template<
unsigned int spacedim>
224 if (ref_data)
delete ref_data;
226 for (
unsigned int sid=0; sid<n_sides_; sid++)
227 delete side_ref_data[sid];
231 template<
unsigned int spacedim>
260 template<
unsigned int spacedim>
264 data.side = cell_side;
288 template<
unsigned int spacedim>
289 template<
unsigned int dim>
301 if (cell().is_valid())
327 for (
unsigned int i=0; i<this->n_points_; i++)
338 for (
unsigned int i=0; i<this->n_points_; i++)
339 data.determinants[i] = det;
343 for (
unsigned int i=0; i<this->n_points_; i++)
344 data.JxW_values[i] = det*this->ref_data->weights[i];
350 arma::mat::fixed<dim,spacedim> ijac;
360 for (
unsigned int i=0; i<this->n_points_; i++)
368 for (
unsigned int i=0; i<this->n_points_; i++)
374 template<
unsigned int spacedim>
375 template<
unsigned int dim>
378 const unsigned int side_idx = side().side_idx();
383 arma::vec::fixed<spacedim> n_cell;
385 n_cell = n_cell/norm(n_cell,2);
386 for (
unsigned int i=0; i<this->n_points_; i++)
395 for (
unsigned int i=0; i<this->n_points_; i++)
396 data.points.set(i) =
Armor::vec<spacedim>( coords*this->side_ref_data[side_idx]->bar_coords[i] );
408 arma::mat::fixed<spacedim,dim> side_coords;
409 arma::mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac;
412 for (
unsigned int n=0; n<dim; n++)
413 for (
unsigned int c=0; c<spacedim; c++)
414 side_coords(c,n) = (*data.side.node(n))[c];
420 for (
unsigned int i=0; i<this->n_points_; i++)
421 data.side_JxW_values[i] = side_det*this->side_ref_data[side_idx]->weights[i];
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
ArmaVec< Type, nr > vec(uint mat_index) const
unsigned int size() const
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
ElementAccessor< spacedim > cell
Iterator to last updated cell.
ElementData(unsigned int size, UpdateFlags flags, unsigned int dim)
Resize the data arrays.
void print()
Print calculated data.
std::vector< double > JxW_values
Transformed quadrature weights.
Armor::array normal_vectors
Normal vectors (spacedim) to the element at the quadrature points lying on a side.
std::vector< double > side_JxW_values
JxW values for sides.
Side side
Iterator to last updated cell side.
Class for computation of data on cell and side.
ElementValues(Quadrature &_quadrature, UpdateFlags _flags, unsigned int dim)
Constructor.
void fill_side_data()
Calculates the mapping data on a side of a cell.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
void fill_data()
Compute data from reference cell and using MappingP1.
Affine mapping between reference and actual cell.
static arma::mat::fixed< spacedim, dim > jacobian(const ElementMap &coords)
arma::mat::fixed< spacedim, dim+1 > ElementMap
static ElementMap element_map(ElementAccessor< 3 > elm)
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Base class for quadrature rules on simplices in arbitrary dimensions.
Quadrature make_from_side(unsigned int sid) const
unsigned int size() const
Returns number of quadrature points.
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
double weight(unsigned int i) const
Returns the ith weight.
Structure for storing the precomputed element data.
RefElementData(unsigned int np)
Resize vectors to size np.
std::vector< arma::vec > bar_coords
Barycentric coordinates of quadrature points.
std::vector< double > weights
Quadrature weights.
void ref_initialize(Quadrature &_quadrature, unsigned int dim)
Initialize ref_data or side_ref_data.
RefElementData * init_ref_data(const Quadrature &q)
Precompute data on reference element.
~RefElementValues()
Correct deallocation of objects created by 'initialize' methods.
RefElementValues(Quadrature &_quadrature, unsigned int dim)
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
static LocalPoint normal_vector(unsigned int sid)
unsigned int elem_idx() const
Returns index of element in Mesh::element_vec_.
unsigned int side_idx() const
Returns local index of the side on the element.
bool is_valid() const
Returns true if the side has assigned element.
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
ArmaMat< double, N, M > mat
void printf(BasicWriter< Char > &w, BasicCStringRef< Char > format, ArgList args)
Basic definitions of numerical quadrature rules.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
@ update_volume_elements
Determinant of the Jacobian.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
@ update_side_JxW_values
Transformed quadrature weight for cell sides.
@ update_inverse_jacobians
Volume element.
@ update_quadrature_points
Transformed quadrature points.