37 template<
int rank,
int elemdim,
int spacedim,
class Value>
41 inline static typename Value::return_type
fe_value(
FEValues<elemdim,3> &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 elemdim,
int spacedim,
class Value>
56 inline static typename Value::return_type
fe_value(
FEValues<elemdim,3> &fe_val,
unsigned int i_dof,
unsigned int i_qp,
unsigned int comp_index)
64 template<
int elemdim,
int spacedim,
class Value>
67 inline static typename Value::return_type
fe_value(
FEValues<elemdim,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 elemdim,
int spacedim,
class Value>
78 inline static typename Value::return_type
fe_value(
FEValues<elemdim,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>
93 template <
int elemdim,
int spacedim,
class Value>
97 WarningOut() <<
"Multiple initialization of FEValueHandler!";
110 template <
int elemdim,
int spacedim,
class Value>
inline 114 point_list.push_back(p);
123 template <
int elemdim,
int spacedim,
class Value>
128 ASSERT_EQ( point_list.size(), value_list.size() ).error();
135 for (
unsigned int k=0; k<point_list.size(); k++) {
142 Value envelope(value_list[k]);
144 for (
unsigned int i=0; i<
dh_->ds()->fe<elemdim>(elm)->n_dofs(); i++) {
152 template <
int elemdim,
int spacedim,
class Value>
156 static const double weight_coefs[] = { 1., 1., 2., 6. };
158 QGauss qgauss(elemdim, order);
161 for(
unsigned i=0; i<qgauss.
size(); ++i) {
162 q_weights[i] = qgauss.
weight(i)*weight_coefs[elemdim];
166 return qgauss.
size();
170 template <
int elemdim,
int spacedim,
class Value>
173 unsigned int ndofs = this->
value_.n_rows() * this->
value_.n_cols();
174 for (
unsigned int k=0; k<ndofs; k++) {
175 indices[k] = (*boundary_dofs_)[ndofs*cell.
idx()+k];
181 template <
int spacedim,
class Value>
185 WarningOut() <<
"Multiple initialization of FEValueHandler!";
194 template <
int spacedim,
class Value>
198 ASSERT_EQ( point_list.size(), value_list.size() ).error();
204 for (
unsigned int k=0; k<point_list.size(); k++) {
205 Value envelope(value_list[k]);
207 for (
unsigned int i=0; i<
dh_->ds()->fe<0>(elm)->n_dofs(); i++) {
214 template <
int spacedim,
class Value>
217 unsigned int ndofs = this->
value_.n_rows() * this->
value_.n_cols();
218 for (
unsigned int k=0; k<ndofs; k++) {
219 indices[k] = (*boundary_dofs_)[ndofs*cell.
idx()+k];
225 template <
int elemdim,
int spacedim,
class Value>
233 #define INSTANCE_VALUE_HANDLER_ALL(dim, spacedim) \ 234 template class FEValueHandler<dim, spacedim, FieldValue<0>::Enum >; \ 235 template class FEValueHandler<dim, spacedim, FieldValue<0>::Integer >; \ 236 template class FEValueHandler<dim, spacedim, FieldValue<0>::Scalar >; \ 237 template class FEValueHandler<dim, spacedim, FieldValue<spacedim>::VectorFixed >; \ 238 template class FEValueHandler<dim, spacedim, FieldValue<spacedim>::TensorFixed >; \ 239 template class FEShapeHandler<0, dim, spacedim, FieldValue<0>::Enum >; \ 240 template class FEShapeHandler<0, dim, spacedim, FieldValue<0>::Integer >; \ 241 template class FEShapeHandler<0, dim, spacedim, FieldValue<0>::Scalar >; \ 242 template class FEShapeHandler<1, dim, spacedim, FieldValue<spacedim>::VectorFixed >; \ 243 template class FEShapeHandler<2, dim, spacedim, FieldValue<spacedim>::TensorFixed >; 245 #define INSTANCE_VALUE_HANDLER(dim) \ 246 INSTANCE_VALUE_HANDLER_ALL(dim,3) unsigned int comp_index
index of component (of vector_value/tensor_value)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
unsigned int comp_index_
Index of component (of vector_value/tensor_value)
RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map) const
void initialize(FEValueInitData init_data)
Initialize data members.
BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map) const
Space< spacedim >::Point Point
unsigned int get_loc_dof_indices(std::vector< LongIdx > &indices) const
Returns the indices of dofs associated to the cell on the local process.
static Value::return_type fe_value(FEValues< elemdim, 3 > &fe_val, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
Cell accessor allow iterate over DOF handler cells.
Class FEValues calculates finite element data on the actual cells such as shape function values...
const FEValuesViews::Tensor< dim, spacedim > & tensor_view(unsigned int i) const
Accessor to tensor values of multicomponent FE.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
double weight(const unsigned int i) const
Returns the ith weight.
Base class for quadrature rules on simplices in arbitrary dimensions.
Symmetric Gauss-Legendre quadrature formulae on simplices.
MappingP1< elemdim, 3 > * map_
Mapping object.
Armor::vec< point_dim > point(const unsigned int i) const
Returns the ith quadrature point.
static Value::return_type fe_value(FEValues< elemdim, 3 > &fe_val, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
void value_list(const std::vector< Point > &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.
Value::return_type r_value_
const unsigned int size() const
Returns number of quadrature points.
Basic definitions of numerical quadrature rules.
std::vector< LongIdx > dof_indices
Array of indexes to data_vec_, used for get/set values.
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
unsigned int ndofs
number of dofs
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.
const FEValuesViews::Vector< dim, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
std::shared_ptr< DOFHandlerMultiDim > dh
DOF handler object.
unsigned int get_dof_indices(const ElementAccessor< 3 > &cell, std::vector< LongIdx > &indices) const
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
MappingP1< elemdim, 3 > * get_mapping()
Return mapping object.
Value value_
Last value, prevents passing large values (vectors) by value.
#define INSTANCE_VALUE_HANDLER(dim)
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
std::shared_ptr< std::vector< LongIdx > > boundary_dofs_
double shape_value(const unsigned int function_no, const unsigned int point_no) override
Return the value of the function_no-th shape function at the point_no-th quadrature point...
VectorMPI data_vec_
Store data of Field.
~FEValueHandler()
Destructor.
static Value::return_type fe_value(FEValues< elemdim, 3 > &fe_val, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining 'warning' record of log.
Calculates finite element data on the actual cell.
static Value::return_type fe_value(FEValues< elemdim, 3 > &fe_val, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
VectorMPI data_vec
Store data of Field.
FEValueHandler()
Constructor.
ElementMap element_map(ElementAccessor< 3 > elm) const
unsigned int n_comp
number of components
Initialization structure of FEValueHandler class.
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)