Flow123d
release_3.0.0-968-gc87a28e79
|
Go to the documentation of this file.
30 template<
unsigned int dim,
unsigned int spacedim>
35 template<
unsigned int dim,
unsigned int spacedim>
51 for (
unsigned int i=0; i<dim; i++)
62 for (
unsigned int i=0; i<q.
size(); i++)
71 template<
unsigned int dim,
unsigned int spacedim>
88 template<
unsigned int dim,
unsigned int spacedim>
95 arma::mat::fixed<spacedim,dim> jac;
103 coords = element_map(cell);
116 for (
unsigned int i=0; i<q.
size(); i++)
126 for (
unsigned int i=0; i<q.
size(); i++)
131 for (
unsigned int i=0; i<q.
size(); i++)
138 arma::mat::fixed<dim,spacedim> ijac;
147 for (
unsigned int i=0; i<q.
size(); i++)
156 for (
unsigned int i=0; i<q.
size(); i++)
161 template<
unsigned int dim,
unsigned int spacedim>
176 coords = element_map(cell);
185 arma::mat::fixed<spacedim,dim> jac = coords*grad;
189 for (
unsigned int i=0; i<q.
size(); i++)
196 for (
unsigned int i=0; i<q.
size(); i++)
203 arma::mat::fixed<dim,spacedim> ijac;
213 for (
unsigned int i=0; i<q.
size(); i++)
219 arma::vec::fixed<spacedim> n_cell;
221 n_cell = n_cell/norm(n_cell,2);
222 for (
unsigned int i=0; i<q.
size(); i++)
232 for (
unsigned int i=0; i<q.
size(); i++)
245 arma::mat::fixed<spacedim,dim> side_coords;
246 arma::mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac;
250 for (
unsigned int n=0; n<dim; n++)
251 for (
unsigned int c=0; c<spacedim; c++)
252 side_coords(c,n) = cell.
side(sid)->
node(n)->point()[c];
253 side_jac = side_coords * grad.submat(0,0,dim-1,dim-2);
258 for (
unsigned int i=0; i<q.
size(); i++)
264 template<
unsigned int dim,
unsigned int spacedim>
268 for (
unsigned int i=0; i<dim+1; i++)
269 coords.col(i) = elm.node(i)->point();
274 template<
unsigned int dim,
unsigned int spacedim>
277 arma::mat::fixed<3, dim> A =
map.cols(1,dim);
278 for(
unsigned int i=0; i < dim; i++ ) {
279 A.col(i) -=
map.col(0);
282 arma::mat::fixed<dim, dim> AtA = A.t()*A;
283 arma::vec::fixed<dim> Atb = A.t()*(point -
map.col(0));
284 arma::vec::fixed<dim+1> bary_coord;
285 bary_coord.rows(1, dim) = arma::solve(AtA, Atb);
286 bary_coord( 0 ) = 1.0 - arma::sum( bary_coord.rows(1,dim) );
290 template<
unsigned int dim,
unsigned int spacedim>
296 template<
unsigned int dim,
unsigned int spacedim>
301 template <
unsigned int dim,
unsigned int spacedim>
304 arma::vec projection = this->project_real_to_unit(point, this->element_map(elm));
static LocalPoint normal_vector(unsigned int sid)
arma::mat::fixed< spacedim, dim+1 > ElementMap
@ update_volume_elements
Determinant of the Jacobian.
const unsigned int size() const
Returns number of quadrature points.
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,...
Mapping data that can be precomputed on the actual cell.
double weight(const unsigned int i) const
Returns the ith weight.
@ update_quadrature_points
Transformed quadrature points.
void fill_fe_values(const ElementAccessor< 3 > &cell, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Calculates the mapping data on the actual cell.
@ update_normal_vectors
Normal vectors.
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
void fill_fe_side_values(const ElementAccessor< 3 > &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.
@ update_inverse_jacobians
Volume element.
static const double epsilon
stabilization parameter
arma::vec::fixed< dim+1 > BaryPoint
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map) const
static BaryPoint clip(const BaryPoint &barycentric)
ElementMap element_map(ElementAccessor< 3 > elm) const
bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
arma::vec::fixed< spacedim > RealPoint
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
BaryPoint clip_to_element(BaryPoint &barycentric)
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
std::vector< double > JxW_values
Transformed quadrature weights.
Basic definitions of numerical quadrature rules.
RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map) const
@ update_JxW_values
Transformed quadrature weights.
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.
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Normal vectors to the element at the quadrature points lying on a side.
std::vector< arma::vec > bar_coords
Auxiliary array of barycentric coordinates of quadrature points.
@ update_jacobians
Volume element.
@ update_side_JxW_values
Transformed quadrature weight for cell sides.
#define ASSERT_LE_DBG(a, b)
Definition of comparative assert macro (Less or Equal) only for debug mode.
NodeAccessor< 3 > node(unsigned int i) const
Base class for quadrature rules on simplices in arbitrary dimensions.
SideIter side(const unsigned int loc_index)
MappingInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)
Initializes the structures and computes static data.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.