Flow123d
JS_before_hm-1727-ga1133b990
|
Go to the documentation of this file.
40 template<
unsigned int spacedim>
50 template<
unsigned int spacedim>
53 unsigned int first_component_idx,
57 for (
unsigned int ip=0; ip<
n_points; ip++)
58 for (
unsigned int id=0;
id<dof_indices.size();
id++)
67 template<
unsigned int spacedim>
68 template<
unsigned int DIM>
88 ASSERT_DBG(fe_sys !=
nullptr).error(
"Mixed system must be represented by FESystem.");
93 unsigned int comp_offset = 0;
94 for (
auto fe : fe_sys->
fe())
122 template<
unsigned int spacedim>
130 template<
unsigned int spacedim>
135 template<
unsigned int spacedim>
136 template<
unsigned int DIM>
142 if (DIM == 0)
return;
144 allocate( q.
size(), _fe, _flags);
145 elm_values = std::make_shared<ElementValues<spacedim> >(q, update_flags, DIM);
151 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
153 fe_values_vec.resize(fe->
fe().size());
154 for (
unsigned int f=0; f<fe->
fe().size(); f++)
155 fe_values_vec[f].initialize(q, *fe->
fe()[f], update_flags);
159 if ( q.dim() == DIM )
161 fe_data = init_fe_data(_fe, q);
163 else if ( q.dim() + 1 == DIM )
166 for (
unsigned int sid = 0; sid < RefElement<DIM>::n_sides; sid++)
172 for (
unsigned int pid = 0; pid < RefElement<DIM>::n_side_permutations; pid++)
173 side_fe_data[sid][pid] = init_fe_data(_fe, q.make_from_side<DIM>(sid,pid));
177 ASSERT_DBG(
false)(q.dim())(DIM).error(
"Dimension mismatch in FEValues::initialize().");
182 template<
unsigned int spacedim>
183 template<
unsigned int DIM>
185 unsigned int n_points,
194 ASSERT_DBG(_fe.
n_components() == spacedim*spacedim).error(
"FETensor must have spacedim*spacedim components.");
197 fe_sys_dofs_.clear();
198 fe_sys_n_components_.clear();
199 fe_sys_n_space_components_.clear();
202 n_points_ = n_points;
205 fe_type_ = _fe.
type_;
207 if (fe_sys !=
nullptr)
209 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
211 fe_sys_dofs_.push_back(fe_sys->
fe_dofs(f));
212 fe_sys_n_components_.push_back(fe_sys->
fe()[f]->n_components());
213 fe_sys_n_space_components_.push_back(fe_sys->
fe()[f]->n_space_components(spacedim));
221 shape_values.resize(n_points_,
vector<double>(n_dofs_*n_components_));
224 shape_gradients.resize(n_points_,
vector<arma::vec::fixed<spacedim> >(n_dofs_*n_components_));
226 views_cache_.initialize(*
this, _fe);
231 template<
unsigned int spacedim>
232 template<
unsigned int DIM>
237 std::shared_ptr<FEInternalData> data = std::make_shared<FEInternalData>(q.
size(), n_dofs_);
240 for (
unsigned int i=0; i<q.
size(); i++)
242 for (
unsigned int j=0; j<n_dofs_; j++)
247 data->ref_shape_values[i][j] = trans(shape_values.row(j));
252 for (
unsigned int i=0; i<q.
size(); i++)
254 for (
unsigned int j=0; j<n_dofs_; j++)
260 data->ref_shape_grads[i][j] = grad;
280 template<
unsigned int spacedim>
282 const unsigned int point_no,
283 const unsigned int comp)
const
288 return shape_gradients[point_no][function_no*n_components_+comp];
497 template<
unsigned int spacedim>
502 this->fill_data_specialized<MapScalar<spacedim>>(elm_values, fe_data);
505 this->fill_data_specialized<MapVector<spacedim>>(elm_values, fe_data);
508 this->fill_data_specialized<MapContravariant<spacedim>>(elm_values, fe_data);
511 this->fill_data_specialized<MapPiola<spacedim>>(elm_values, fe_data);
514 this->fill_data_specialized<MapTensor<spacedim>>(elm_values, fe_data);
517 this->fill_data_specialized<MapSystem<spacedim>>(elm_values, fe_data);
520 ASSERT(
false).error(
"Not implemented.");
526 template<
unsigned int spacedim>
527 template<
class MapType>
530 map_type.fill_values_vec(*
this, elm_values, fe_data);
532 map_type.update_values(*
this, elm_values, fe_data);
534 map_type.update_gradients(*
this, elm_values, fe_data);
540 template<
unsigned int spacedim>
545 if (!elm_values->cell().is_valid() ||
546 elm_values->cell() != cell)
548 elm_values->reinit(cell);
551 fill_data(*elm_values, *fe_data);
555 template<
unsigned int spacedim>
560 if (!elm_values->side().is_valid() ||
561 elm_values->side() != cell_side)
563 elm_values->reinit(cell_side);
567 const unsigned int pid = elm_values->side().element()->permutation_idx(sid);
570 fill_data(*elm_values, *side_fe_data[sid][pid]);
581 fv[0].initialize(quadrature[0], *fe[0_d], flags);
582 fv[1].initialize(quadrature[1], *fe[1_d], flags);
583 fv[2].initialize(quadrature[2], *fe[2_d], flags);
584 fv[3].initialize(quadrature[3], *fe[3_d], flags);
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
~FEValues()
Correct deallocation of objects created by 'initialize' methods.
@ update_gradients
Shape function gradients.
Helper class allows update values and gradients of FEValues of FEVectorPiola type.
unsigned int n_points() const
Returns the number of quadrature points.
arma::vec::fixed< dim > shape_grad(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the comp-th component of the gradient of the i-th shape function at the point p on the ref...
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
FEValues()
Default constructor with postponed initialization.
static 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,...
@ update_values
Shape function values.
unsigned int n_dofs_
Number of finite element dofs.
unsigned int size() const
Returns number of quadrature points.
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
void initialize(const FEValues &fv, const FiniteElement< DIM > &fe)
Class FESystem for compound finite elements.
Helper class allows update values and gradients of FEValues of FEVector type.
Calculates finite element data on the actual cell.
unsigned int dim_
Dimension of reference space.
unsigned int n_components() const
Returns numer of components of the basis function.
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
void fill_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure.
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
virtual UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
Structure for storing the precomputed finite element data.
ArmaMat< double, N, M > mat
FEType type_
Type of FiniteElement.
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Abstract class for description of finite elements.
Helper class allows update values and gradients of FEValues of FETensor type.
std::array< QGauss, 4 > array
unsigned int side_idx() const
Returns local index of the side on the element.
Helper class allows update values and gradients of FEValues of FEMixedSystem type.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
Class for computation of data on cell and side.
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
unsigned int n_points_
Number of integration points.
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Helper class allows update values and gradients of FEValues of FEScalar type.
std::vector< unsigned int > fe_dofs(unsigned int fe_index)
Return dof indices belonging to given sub-FE.
void allocate(unsigned int n_points, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
Basic definitions of numerical quadrature rules.
void fill_data_specialized(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure....
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
arma::vec::fixed< spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
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.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Base class for quadrature rules on simplices in arbitrary dimensions.
Compound finite element on dim dimensional simplex.
std::shared_ptr< FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
Helper class allows update values and gradients of FEValues of FEVectorContravariant type.
unsigned int n_space_components(unsigned int spacedim)
Number of components of FE in a mapped space with dimension spacedim.
FEInternalData(unsigned int np, unsigned int nd)
unsigned int n_points
Number of quadrature points.
double shape_value(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the value of the comp-th component of the i-th shape function at the point p on the refere...