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)
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>
111 static const double weight_coefs[] = { 1., 1., 2., 6. };
113 QGauss qgauss(elemdim, order);
116 for(
unsigned i=0; i<qgauss.
size(); ++i) {
117 q_weights[i] = qgauss.
weight(i)*weight_coefs[elemdim];
121 return qgauss.
size();
125 template <
int spacedim,
class Value>
129 WarningOut() <<
"Multiple initialization of FEValueHandler!";
139 template <
int spacedim,
class Value>
148 for (
unsigned int k=0; k<point_list.
size(); k++) {
149 Value envelope(value_list[k]);
152 envelope(i / envelope.n_cols(), i % envelope.n_rows()) +=
data_vec_.
get(loc_dofs[i]);
158 template <
int elemdim,
int spacedim,
class Value>
165 #define INSTANCE_VALUE_HANDLER_ALL(dim, spacedim) \
166 template class FEValueHandler<dim, spacedim, FieldValue<0>::Enum >; \
167 template class FEValueHandler<dim, spacedim, FieldValue<0>::Integer >; \
168 template class FEValueHandler<dim, spacedim, FieldValue<0>::Scalar >; \
169 template class FEValueHandler<dim, spacedim, FieldValue<spacedim>::VectorFixed >; \
170 template class FEValueHandler<dim, spacedim, FieldValue<spacedim>::TensorFixed >;
172 #define INSTANCE_VALUE_HANDLER(dim) \
173 INSTANCE_VALUE_HANDLER_ALL(dim,3)
~FEValueHandler()
Destructor.
arma::Col< IntIdx > LocDofVec
unsigned int range_end
Holds end of component range of evaluation.
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.
Class FEValues calculates finite element data on the actual cells such as 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.
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definitions of particular quadrature rules on simplices.
unsigned int range_end_
End of dof range of actual component.
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) only for debug mode.
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.
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.
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
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
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
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.
unsigned int range_begin
Holds begin of component range of evaluation.
#define INSTANCE_VALUE_HANDLER(dim)
unsigned int n_comp
number of components