Flow123d
JS_before_hm-1804-gf2ad740aa
|
Go to the documentation of this file.
29 template<
unsigned int dim,
unsigned int spacedim>
47 template<
unsigned int dim,
unsigned int spacedim>
51 for (
unsigned int i=0; i<dim+1; i++)
52 coords.col(i) = *elm.node(i);
57 template<
unsigned int dim,
unsigned int spacedim>
60 arma::mat::fixed<spacedim,dim> jac;
61 for (
unsigned int i=0; i<spacedim; i++)
62 for (
unsigned int j=0; j<dim; j++)
63 jac(i,j) = coords(i,j+1) - coords(i,0);
68 template<
unsigned int dim,
unsigned int spacedim>
71 arma::mat::fixed<3, dim> A =
map.cols(1,dim);
72 for(
unsigned int i=0; i < dim; i++ ) {
73 A.col(i) -=
map.col(0);
76 arma::mat::fixed<dim, dim> AtA = A.t()*A;
77 arma::vec::fixed<dim> Atb = A.t()*(point -
map.col(0));
78 arma::vec::fixed<dim+1> bary_coord;
79 bary_coord.rows(1, dim) = arma::solve(AtA, Atb);
80 bary_coord( 0 ) = 1.0 - arma::sum( bary_coord.rows(1,dim) );
84 template<
unsigned int dim,
unsigned int spacedim>
90 template<
unsigned int dim,
unsigned int spacedim>
95 template <
unsigned int dim,
unsigned int spacedim>
98 arma::vec projection = project_real_to_unit(point, element_map(elm));
arma::mat::fixed< spacedim, dim+1 > ElementMap
@ update_volume_elements
Determinant of the Jacobian.
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map)
@ update_normal_vectors
Normal vectors.
@ update_inverse_jacobians
Volume element.
static const double epsilon
stabilization parameter
static RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map)
arma::vec::fixed< dim+1 > BaryPoint
static BaryPoint clip(const BaryPoint &barycentric)
static bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
arma::vec::fixed< spacedim > RealPoint
static arma::mat::fixed< spacedim, dim > jacobian(const ElementMap &coords)
static ElementMap element_map(ElementAccessor< 3 > elm)
static BaryPoint clip_to_element(BaryPoint &barycentric)
Affine mapping between reference and actual cell.
Basic definitions of numerical quadrature rules.
@ 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.
@ update_jacobians
Volume element.