Flow123d
JS_before_hm-1820-g636551b33
|
Go to the documentation of this file.
37 template<
int rank,
int spacedim,
class Value>
41 inline static typename Value::return_type
fe_value(
FEValues<spacedim> &fe_val,
unsigned int i_dof,
unsigned int i_qp,
unsigned int comp_index)
43 ASSERT(
false).error(
"Unsupported format of FieldFE!\n");
44 typename Value::return_type ret;
53 template<
int spacedim,
class Value>
58 return fe_val.
scalar_view(comp_index).value(i_dof, i_qp);
64 template<
int spacedim,
class Value>
67 inline static typename Value::return_type
fe_value(
FEValues<3> &fe_val,
unsigned int i_dof,
unsigned int i_qp,
unsigned int comp_index)
69 return fe_val.
vector_view(comp_index).value(i_dof, i_qp);
75 template<
int spacedim,
class Value>
78 inline static typename Value::return_type
fe_value(
FEValues<3> &fe_val,
unsigned int i_dof,
unsigned int i_qp,
unsigned int comp_index)
80 return fe_val.
tensor_view(comp_index).value(i_dof, i_qp);
86 template <
int elemdim,
int spacedim,
class Value>
92 template <
int elemdim,
int spacedim,
class Value>
96 WarningOut() <<
"Multiple initialization of FEValueHandler!";
100 value_.set_n_comp(init_data.
n_comp);
107 template <
int elemdim,
int spacedim,
class Value>
inline
113 v_list.push_back(r_value_);
114 this->value_list(point_list, elm, v_list);
115 this->r_value_ = v_list[0];
116 return this->r_value_;
120 template <
int elemdim,
int spacedim,
class Value>
128 if (boundary_dofs_) loc_dofs = this->get_loc_dof_indices(elm.
idx());
134 for (
unsigned int k=0; k<point_list.
size(); k++)
140 for (
unsigned int k=0; k<point_list.
size(); k++) {
141 Value envelope(value_list[k]);
143 for (
unsigned int i=this->range_begin_, i_dof=0; i<this->range_end_; i++, i_dof++) {
144 value_list[k] += data_vec_.get(loc_dofs[i])
151 template <
int elemdim,
int spacedim,
class Value>
155 static const double weight_coefs[] = { 1., 1., 2., 6. };
157 QGauss qgauss(elemdim, order);
160 for(
unsigned i=0; i<qgauss.
size(); ++i) {
161 q_weights[i] = qgauss.
weight(i)*weight_coefs[elemdim];
165 return qgauss.
size();
169 template <
int spacedim,
class Value>
173 WarningOut() <<
"Multiple initialization of FEValueHandler!";
183 template <
int spacedim,
class Value>
194 for (
unsigned int k=0; k<point_list.
size(); k++) {
198 envelope(i / envelope.n_cols(), i % envelope.n_rows()) +=
data_vec_.
get(loc_dofs[i]);
204 template <
int elemdim,
int spacedim,
class Value>
211 #define INSTANCE_VALUE_HANDLER_ALL(dim, spacedim) \
212 template class FEValueHandler<dim, spacedim, FieldValue<0>::Enum >; \
213 template class FEValueHandler<dim, spacedim, FieldValue<0>::Integer >; \
214 template class FEValueHandler<dim, spacedim, FieldValue<0>::Scalar >; \
215 template class FEValueHandler<dim, spacedim, FieldValue<spacedim>::VectorFixed >; \
216 template class FEValueHandler<dim, spacedim, FieldValue<spacedim>::TensorFixed >;
218 #define INSTANCE_VALUE_HANDLER(dim) \
219 INSTANCE_VALUE_HANDLER_ALL(dim,3)
~FEValueHandler()
Destructor.
arma::Col< IntIdx > LocDofVec
unsigned int range_end
Holds end of component range of evaluation.
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
LocDofVec get_loc_dof_indices(unsigned int cell_idx) const
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh....
Space< spacedim >::Point Point
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
void initialize(FEValueInitData init_data)
Initialize data members.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
ArmaVec< Type, nr > vec(uint mat_index) const
Class FEValues calculates finite element data on the actual cells such as shape function values,...
@ update_values
Shape function values.
unsigned int size() const
Returns number of quadrature points.
Calculates finite element data on the actual cell.
MixedPtr< FiniteElement > mixed_fe
FiniteElement objects of all dimensions.
static RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map)
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
ArmaMat< double, N, M > mat
std::shared_ptr< DOFHandlerMultiDim > dh
DOF handler object.
Definitions of particular quadrature rules on simplices.
unsigned int range_end_
End of dof range of actual component.
static LocalPoint bary_to_local(const BaryPoint &bp)
Converts from barycentric to local coordinates.
Armor::Array< double >::ArrayMatSet set(uint i)
double weight(unsigned int i) const
Returns the ith weight.
static Value::return_type fe_value(FEValues< 3 > &fe_val, unsigned int i_dof, unsigned int i_qp, FMT_UNUSED unsigned int comp_index)
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
static Value::return_type fe_value(FEValues< 3 > &fe_val, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
VectorMPI data_vec_
Store data of Field.
static Value::return_type fe_value(FEValues< 3 > &fe_val, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
Value value_
Last value, prevents passing large values (vectors) by value.
unsigned int compute_quadrature(std::vector< arma::vec::fixed< 3 >> &q_points, std::vector< double > &q_weights, const ElementAccessor< spacedim > &elm, unsigned int order=3)
Compute real coordinates and weights (use QGauss) for given element.
ArrayMatSet set(uint index)
unsigned int size() const
Initialization structure of FEValueHandler class.
const FEValuesViews::Tensor< spacedim > & tensor_view(unsigned int i) const
Accessor to tensor values of multicomponent FE.
Symmetric Gauss-Legendre quadrature formulae on simplices.
FEValueHandler()
Constructor.
static ElementMap element_map(ElementAccessor< 3 > elm)
Cell accessor allow iterate over DOF handler cells.
#define WarningOut()
Macro defining 'warning' record of log.
const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Affine mapping between reference and actual cell.
const FEValuesViews::Scalar< spacedim > & scalar_view(unsigned int i) const
Accessor to scalar values of multicomponent FE.
Basic definitions of numerical quadrature rules.
unsigned int range_begin_
Begin of dof range of actual component.
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
VectorMPI data_vec
Store data of Field.
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
static Value::return_type fe_value(FEValues< spacedim > &fe_val, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
double get(unsigned int pos) const
Return value on given position.
void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Returns std::vector of scalar values in several points at once.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
unsigned int range_begin
Holds begin of component range of evaluation.
#define INSTANCE_VALUE_HANDLER(dim)
Base class for quadrature rules on simplices in arbitrary dimensions.
typename arma::Col< Type >::template fixed< nr > ArmaVec
std::shared_ptr< std::vector< IntIdx > > boundary_dofs_
unsigned int n_comp
number of components