Flow123d
JS_before_hm-1713-g943b6bc16
|
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),
162 n_side_permutations_(_quadrature.dim() +1 == dim ? (dim+1)*(2*dim*dim-5*dim+6)/6 : 0),
165 data(n_points_, update_each(_flags), dim)
167 if (dim == 0)
return;
168 if ( _quadrature.
dim() == dim )
173 else if ( _quadrature.
dim() + 1 == dim )
176 for (
unsigned int sid = 0; sid <
n_sides_; sid++)
194 ASSERT(
false)(dim).error(
"Unsupported dimension.\n");
204 template<
unsigned int spacedim>
207 if (ref_data)
delete ref_data;
209 for (
unsigned int sid=0; sid<n_sides_; sid++)
210 for (
unsigned int pid=0; pid<n_side_permutations_; pid++)
211 delete side_ref_data[sid][pid];
215 template<
unsigned int spacedim>
234 ASSERT(
false)(dim_).error(
"Unsupported dimension.\n");
240 template<
unsigned int spacedim>
244 data.side = cell_side;
262 ASSERT(
false)(dim_).error(
"Unsupported dimension.\n");
268 template<
unsigned int spacedim>
269 template<
unsigned int dim>
281 if (cell().is_valid())
298 for (
unsigned int i=0; i<n_points_; i++)
309 for (
unsigned int i=0; i<n_points_; i++)
310 data.determinants[i] = det;
314 for (
unsigned int i=0; i<n_points_; i++)
315 data.JxW_values[i] = det*ref_data->weights[i];
321 arma::mat::fixed<dim,spacedim> ijac;
330 for (
unsigned int i=0; i<n_points_; i++)
338 for (
unsigned int i=0; i<n_points_; i++)
344 template<
unsigned int spacedim>
345 template<
unsigned int dim>
348 const unsigned int side_idx = side().side_idx();
349 const unsigned int perm_idx = side().element()->permutation_idx(side_idx);
354 arma::vec::fixed<spacedim> n_cell;
356 n_cell = n_cell/norm(n_cell,2);
357 for (
unsigned int i=0; i<n_points_; i++)
366 for (
unsigned int i=0; i<n_points_; i++)
367 data.points.set(i) =
Armor::vec<spacedim>( coords*side_ref_data[side_idx][perm_idx]->bar_coords[i] );
379 arma::mat::fixed<spacedim,dim> side_coords;
380 arma::mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac;
383 for (
unsigned int n=0; n<dim; n++)
384 for (
unsigned int c=0; c<spacedim; c++)
385 side_coords(c,n) = (*data.side.node(n))[c];
391 for (
unsigned int i=0; i<n_points_; i++)
392 data.side_JxW_values[i] = side_det*side_ref_data[side_idx][perm_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< std::vector< RefElementData * > > side_ref_data
Data on reference element (for each side and its permutation).
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.
const unsigned int n_side_permutations_
Number of permutations of points on side of reference cell.
Quadrature make_from_side(unsigned int sid, unsigned int pid) const
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.