Flow123d  release_2.2.0-914-gf1a3a4f
fe_value_handler.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file fe_value_handler.cc
15  * @brief
16  */
17 
19 #include "fem/mapping_p1.hh"
20 #include "fem/fe_values.hh"
21 #include "quadrature/quadrature.hh"
22 #include "mesh/bounding_box.hh"
23 #include "fem/fe_values_views.hh"
24 
25 
26 /**
27  * Helper class, allow to simplify computing value of FieldFE.
28  *
29  * Use correct method FEValues<...>::shape_xxx given with Value::rank_.
30  * Practical use have only instances with rank template parameters 0 and 1 (Scalar and Vector Fields, see below).
31  */
32 template<int rank, int elemdim, int spacedim, class Value>
34 public:
35 
36  inline static typename Value::return_type fe_value(FEValues<elemdim,3> &fe_val, unsigned int i_dof, unsigned int i_qp)
37  {
38  ASSERT(false).error("Unsupported format of FieldFE!\n");
39  typename Value::return_type ret;
40  Value val(ret);
41  val.zeros();
42  return ret;
43  }
44 };
45 
46 
47 /// Partial template specialization of FEShapeHandler for scalar fields
48 template<int elemdim, int spacedim, class Value>
49 class FEShapeHandler<0, elemdim, spacedim, Value> {
50 public:
51  inline static typename Value::return_type fe_value(FEValues<elemdim,3> &fe_val, unsigned int i_dof, unsigned int i_qp)
52  {
53  return fe_val.shape_value(i_dof, i_qp);
54  }
55 };
56 
57 
58 /// Partial template specialization of FEShapeHandler for vector fields
59 template<int elemdim, int spacedim, class Value>
60 class FEShapeHandler<1, elemdim, spacedim, Value> {
61 public:
62  inline static typename Value::return_type fe_value(FEValues<elemdim,3> &fe_val, unsigned int i_dof, unsigned int i_qp)
63  {
64  return fe_val.vector_view(0).value(i_dof, i_qp);
65  }
66 };
67 
68 
69 
70 template <int elemdim, int spacedim, class Value>
72 : value_(r_value_),
73  map_(nullptr)
74 {}
75 
76 
77 template <int elemdim, int spacedim, class Value>
79 {
80  ASSERT_EQ(dof_indices.size(), 0).error("Multiple initialization.");
81 
82  dh_ = init_data.dh;
83  data_vec_ = init_data.data_vec;
84  dof_indices.resize(init_data.ndofs);
85  value_.set_n_comp(init_data.n_comp);
86 
87  if (map == nullptr) {
88  // temporary solution - these objects will be set through FieldCommon
89  map_ = new MappingP1<elemdim,3>();
90  } else {
91  map_ = map;
92  }
93 }
94 
95 
96 template <int elemdim, int spacedim, class Value> inline
97 typename Value::return_type const &FEValueHandler<elemdim, spacedim, Value>::value(const Point &p, const ElementAccessor<spacedim> &elm)
98 {
99  std::vector<Point> point_list;
100  point_list.push_back(p);
102  v_list.push_back(r_value_);
103  this->value_list(point_list, elm, v_list);
104  this->r_value_ = v_list[0];
105  return this->r_value_;
106 }
107 
108 
109 template <int elemdim, int spacedim, class Value>
112 {
113  ASSERT_PTR(map_).error();
114  ASSERT_EQ( point_list.size(), value_list.size() ).error();
115 
116  DOFHandlerBase::CellIterator cell = dh_->mesh()->element( elm.idx() );
117  dh_->get_dof_indices(cell, dof_indices);
118 
119  arma::mat map_mat = map_->element_map(*elm.element());
120  for (unsigned int k=0; k<point_list.size(); k++) {
121  Quadrature<elemdim> quad(1);
122  quad.set_point(0, RefElement<elemdim>::bary_to_local(map_->project_real_to_unit(point_list[k], map_mat)));
123 
124  FEValues<elemdim,3> fe_values(*this->get_mapping(), quad, *dh_->fe<elemdim>(), update_values);
125  fe_values.reinit(cell);
126 
127  Value envelope(value_list[k]);
128  envelope.zeros();
129  for (unsigned int i=0; i<dh_->fe<elemdim>()->n_dofs(); i++) {
130  value_list[k] += (*data_vec_)[dof_indices[i]]
132  }
133  }
134 }
135 
136 
137 template <int elemdim, int spacedim, class Value>
139 {
140  ASSERT_PTR(map_).error();
141 
142  arma::vec projection = map_->project_real_to_unit(point, map_->element_map(elm));
143  return (projection.min() >= -BoundingBox::epsilon);
144 }
145 
146 
147 template <int elemdim, int spacedim, class Value>
149 {}
150 
151 
152 // Instantiations of FEValueHandler and FEShapeHandler
153 #define INSTANCE_VALUE_HANDLER_ALL(dim, spacedim) \
154 template class FEValueHandler<dim, spacedim, FieldValue<0>::Enum >; \
155 template class FEValueHandler<dim, spacedim, FieldValue<0>::Integer >; \
156 template class FEValueHandler<dim, spacedim, FieldValue<0>::Scalar >; \
157 template class FEValueHandler<dim, spacedim, FieldValue<spacedim>::VectorFixed >; \
158 template class FEValueHandler<dim, spacedim, FieldValue<spacedim>::TensorFixed >; \
159 template class FEShapeHandler<0, dim, spacedim, FieldValue<0>::Enum >; \
160 template class FEShapeHandler<0, dim, spacedim, FieldValue<0>::Integer >; \
161 template class FEShapeHandler<0, dim, spacedim, FieldValue<0>::Scalar >; \
162 template class FEShapeHandler<1, dim, spacedim, FieldValue<spacedim>::VectorFixed >; \
163 template class FEShapeHandler<2, dim, spacedim, FieldValue<spacedim>::TensorFixed >;
164 
165 #define INSTANCE_VALUE_HANDLER(dim) \
166 INSTANCE_VALUE_HANDLER_ALL(dim,2) \
167 INSTANCE_VALUE_HANDLER_ALL(dim,3)
168 
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:87
const Element * element() const
Definition: accessors.hh:86
BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map) const
Definition: mapping_p1.cc:275
void initialize(FEValueInitData init_data, MappingP1< elemdim, 3 > *map=nullptr)
Initialize data members.
Space< spacedim >::Point Point
static Value::return_type fe_value(FEValues< elemdim, 3 > &fe_val, unsigned int i_dof, unsigned int i_qp)
Class FEValues calculates finite element data on the actual cells such as shape function values...
static const double epsilon
stabilization parameter
Definition: bounding_box.hh:57
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:32
MappingP1< elemdim, 3 > * map_
Mapping object.
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_
VectorSeqDouble * data_vec
Store data of Field.
std::vector< IdxInt > dof_indices
Array of indexes to data_vec_, used for get/set values.
Basic definitions of numerical quadrature rules.
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
unsigned int ndofs
number of dofs
const FEValuesViews::Vector< dim, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:349
VectorSeqDouble * data_vec_
Store data of Field.
std::shared_ptr< DOFHandlerMultiDim > dh
DOF handler object.
ElementMap element_map(const Element &elm) const
Definition: mapping_p1.cc:265
double shape_value(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:169
static Value::return_type fe_value(FEValues< elemdim, 3 > &fe_val, unsigned int i_dof, unsigned int i_qp)
MappingP1< elemdim, 3 > * get_mapping()
Return mapping object.
Value value_
Last value, prevents passing large values (vectors) by value.
#define INSTANCE_VALUE_HANDLER(dim)
void set_point(const unsigned int i, const arma::vec::fixed< dim > &p)
Sets individual quadrature point coordinates.
Definition: quadrature.hh:155
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:335
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:233
~FEValueHandler()
Destructor.
Calculates finite element data on the actual cell.
Definition: fe_values.hh:450
bool contains_point(arma::vec point, Element &elm)
Test if element contains given point.
FEValueHandler()
Constructor.
unsigned int n_comp
number of components
Initialization structure of FEValueHandler class.
unsigned int idx() const
Definition: accessors.hh:113
static Value::return_type fe_value(FEValues< elemdim, 3 > &fe_val, unsigned int i_dof, unsigned int i_qp)
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:327