Flow123d
JS_before_hm-2198-g122e1f2e2
|
Go to the documentation of this file.
45 template<
unsigned int spacedim>
61 template<
unsigned int spacedim>
64 if (cell.is_valid() || side.is_valid())
67 printf(
"cell %d dim %d ", cell.elm_idx(), cell.dim());
68 else if (side.is_valid())
69 printf(
"cell %d dim %d side %d ", side.elem_idx(), side.dim(), side.side_idx());
72 for (
auto d : determinants)
printf(
"%g ", d);
printf(
"]");
78 for (
auto j : side_JxW_values)
printf(
"%g ", j);
printf(
"]");
81 for (
unsigned int i=0; i<normal_vectors.size(); i++)
83 auto n = normal_vectors.vec<spacedim>(i);
85 for (
unsigned int c=0; c<spacedim; c++)
printf(
"%g ", n[c]);
96 template<
unsigned int spacedim>
103 for (
unsigned int i=0; i<q.
size(); i++)
117 ASSERT(
false)(q.
dim()).error(
"Unsupported dimension.\n");
128 template<
unsigned int spacedim>
146 ASSERT(
false)(dim_).error(
"Unsupported dimension.\n");
154 template<
unsigned int spacedim>
160 n_points_(_quadrature.size()),
161 n_sides_(_quadrature.dim() == dim ? 0 : dim+1),
163 side_ref_data(n_sides_),
164 data(n_points_, update_each(_flags), dim)
166 if (dim == 0)
return;
167 if ( _quadrature.
dim() == dim )
172 else if ( _quadrature.
dim() + 1 == dim )
175 for (
unsigned int sid = 0; sid <
n_sides_; sid++)
191 ASSERT(
false)(dim).error(
"Unsupported dimension.\n");
200 template<
unsigned int spacedim>
203 if (ref_data)
delete ref_data;
205 for (
unsigned int sid=0; sid<n_sides_; sid++)
206 delete side_ref_data[sid];
210 template<
unsigned int spacedim>
229 ASSERT(
false)(dim_).error(
"Unsupported dimension.\n");
235 template<
unsigned int spacedim>
239 data.side = cell_side;
257 ASSERT(
false)(dim_).error(
"Unsupported dimension.\n");
263 template<
unsigned int spacedim>
264 template<
unsigned int dim>
276 if (cell().is_valid())
293 for (
unsigned int i=0; i<n_points_; i++)
304 for (
unsigned int i=0; i<n_points_; i++)
305 data.determinants[i] = det;
309 for (
unsigned int i=0; i<n_points_; i++)
310 data.JxW_values[i] = det*ref_data->weights[i];
316 arma::mat::fixed<dim,spacedim> ijac;
325 for (
unsigned int i=0; i<n_points_; i++)
333 for (
unsigned int i=0; i<n_points_; i++)
339 template<
unsigned int spacedim>
340 template<
unsigned int dim>
343 const unsigned int side_idx = side().side_idx();
348 arma::vec::fixed<spacedim> n_cell;
350 n_cell = n_cell/norm(n_cell,2);
351 for (
unsigned int i=0; i<n_points_; i++)
360 for (
unsigned int i=0; i<n_points_; i++)
373 arma::mat::fixed<spacedim,dim> side_coords;
374 arma::mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac;
377 for (
unsigned int n=0; n<dim; n++)
378 for (
unsigned int c=0; c<spacedim; c++)
379 side_coords(c,n) = (*data.side.node(n))[c];
385 for (
unsigned int i=0; i<n_points_; i++)
386 data.side_JxW_values[i] = side_det*side_ref_data[side_idx]->weights[i];
#define OLD_ASSERT_EQUAL(a, b)
static LocalPoint normal_vector(unsigned int sid)
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
arma::mat::fixed< spacedim, dim+1 > ElementMap
@ update_volume_elements
Determinant of the Jacobian.
std::vector< double > weights
Quadrature weights.
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
RefElementData(unsigned int np)
Resize vectors to size np.
RefElementData * init_ref_data(const Quadrature &q)
Precompute data on reference element.
Structure for storing the precomputed element data.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
std::vector< RefElementData * > side_ref_data
Data on reference element (for each side ).
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
unsigned int size() const
Returns number of quadrature points.
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
RefElementData * ref_data
Data on reference element.
@ update_quadrature_points
Transformed quadrature points.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
@ update_normal_vectors
Normal vectors.
void fill_data()
Compute data from reference cell and using MappingP1.
~ElementValues()
Correct deallocation of objects created by 'initialize' methods.
@ update_inverse_jacobians
Volume element.
ElementValues(Quadrature &_quadrature, UpdateFlags _flags, unsigned int dim)
Constructor.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
ArmaMat< double, N, M > mat
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
void fill_side_data()
Calculates the mapping data on a side of a cell.
double weight(unsigned int i) const
Returns the ith weight.
Class for computation of data on cell and side.
const unsigned int n_sides_
Number of sides in reference cell.
static arma::mat::fixed< spacedim, dim > jacobian(const ElementMap &coords)
static ElementMap element_map(ElementAccessor< 3 > elm)
void print()
Print calculated data.
ElementData(unsigned int size, UpdateFlags flags, unsigned int dim)
Resize the data arrays.
void printf(BasicWriter< Char > &w, BasicCStringRef< Char > format, ArgList args)
Affine mapping between reference and actual cell.
Basic definitions of numerical quadrature rules.
@ update_JxW_values
Transformed quadrature weights.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
@ update_jacobians
Volume element.
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
@ update_side_JxW_values
Transformed quadrature weight for cell sides.
Base class for quadrature rules on simplices in arbitrary dimensions.
std::vector< arma::vec > bar_coords
Barycentric coordinates of quadrature points.
Quadrature make_from_side(unsigned int sid) const