Flow123d  release_3.0.0-688-g6b683cf
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 "la/vector_mpi.hh"
20 #include "fem/mapping_p1.hh"
21 #include "fem/fe_values.hh"
22 #include "quadrature/quadrature.hh"
24 #include "mesh/bounding_box.hh"
25 #include "mesh/accessors.hh"
26 #include "fem/fe_values_views.hh"
27 #include "fem/dh_cell_accessor.hh"
28 
29 
30 /**
31  * Helper class, allow to simplify computing value of FieldFE.
32  *
33  * Use correct method FEValues<...>::shape_xxx given with Value::rank_.
34  * Practical use have only instances with rank template parameters 0 and 1 (Scalar and Vector Fields, see below).
35  */
36 template<int rank, int elemdim, int spacedim, class Value>
38 public:
39 
40  inline static typename Value::return_type fe_value(FEValues<elemdim,3> &fe_val, unsigned int i_dof, unsigned int i_qp)
41  {
42  ASSERT(false).error("Unsupported format of FieldFE!\n");
43  typename Value::return_type ret;
44  Value val(ret);
45  val.zeros();
46  return ret;
47  }
48 };
49 
50 
51 /// Partial template specialization of FEShapeHandler for scalar fields
52 template<int elemdim, int spacedim, class Value>
53 class FEShapeHandler<0, elemdim, spacedim, Value> {
54 public:
55  inline static typename Value::return_type fe_value(FEValues<elemdim,3> &fe_val, unsigned int i_dof, unsigned int i_qp)
56  {
57  return fe_val.shape_value(i_dof, i_qp);
58  }
59 };
60 
61 
62 /// Partial template specialization of FEShapeHandler for vector fields
63 template<int elemdim, int spacedim, class Value>
64 class FEShapeHandler<1, elemdim, spacedim, Value> {
65 public:
66  inline static typename Value::return_type fe_value(FEValues<elemdim,3> &fe_val, unsigned int i_dof, unsigned int i_qp)
67  {
68  return fe_val.vector_view(0).value(i_dof, i_qp);
69  }
70 };
71 
72 
73 /// Partial template specialization of FEShapeHandler for tensor fields
74 template<int elemdim, int spacedim, class Value>
75 class FEShapeHandler<2, elemdim, spacedim, Value> {
76 public:
77  inline static typename Value::return_type fe_value(FEValues<elemdim,3> &fe_val, unsigned int i_dof, unsigned int i_qp)
78  {
79  return fe_val.tensor_view(0).value(i_dof, i_qp);
80  }
81 };
82 
83 
84 
85 template <int elemdim, int spacedim, class Value>
87 : value_(r_value_),
88  map_(nullptr)
89 {}
90 
91 
92 template <int elemdim, int spacedim, class Value>
94 {
95  if (dof_indices.size() > 0)
96  WarningOut() << "Multiple initialization of FEValueHandler!";
97 
98  dh_ = init_data.dh;
99  data_vec_ = init_data.data_vec;
100  dof_indices.resize(init_data.ndofs);
101  value_.set_n_comp(init_data.n_comp);
102 
103  if (map == nullptr) {
104  // temporary solution - these objects will be set through FieldCommon
105  map_ = new MappingP1<elemdim,3>();
106  } else {
107  map_ = map;
108  }
109 }
110 
111 
112 template <int elemdim, int spacedim, class Value> inline
113 typename Value::return_type const &FEValueHandler<elemdim, spacedim, Value>::value(const Point &p, const ElementAccessor<spacedim> &elm)
114 {
115  std::vector<Point> point_list;
116  point_list.push_back(p);
118  v_list.push_back(r_value_);
119  this->value_list(point_list, elm, v_list);
120  this->r_value_ = v_list[0];
121  return this->r_value_;
122 }
123 
124 
125 template <int elemdim, int spacedim, class Value>
128 {
129  ASSERT_PTR(map_).error();
130  ASSERT_EQ( point_list.size(), value_list.size() ).error();
131 
132  ElementAccessor<3> cell = dh_->mesh()->element_accessor( elm.mesh_idx() ); // non-const ElementAccessor
133  if (boundary_dofs_) this->get_dof_indices( elm, dof_indices);
134  else dh_->get_dof_indices( cell, dof_indices );
135 
136  arma::mat map_mat = map_->element_map(elm);
137  for (unsigned int k=0; k<point_list.size(); k++) {
138  Quadrature<elemdim> quad(1);
139  quad.set_point(0, RefElement<elemdim>::bary_to_local(map_->project_real_to_unit(point_list[k], map_mat)));
140 
141  FEValues<elemdim,3> fe_values(*this->get_mapping(), quad, *dh_->fe<elemdim>(cell), update_values);
142  fe_values.reinit( const_cast<ElementAccessor<spacedim> &>(elm) );
143 
144  Value envelope(value_list[k]);
145  envelope.zeros();
146  for (unsigned int i=0; i<dh_->fe<elemdim>(cell)->n_dofs(); i++) {
147  value_list[k] += (*data_vec_)[dof_indices[i]]
149  }
150  }
151 }
152 
153 
154 template <int elemdim, int spacedim, class Value>
155 unsigned int FEValueHandler<elemdim, spacedim, Value>::compute_quadrature(std::vector<arma::vec::fixed<3>> & q_points, std::vector<double> & q_weights,
156  const ElementAccessor<spacedim> &ele, unsigned int order)
157 {
158  static const double weight_coefs[] = { 1., 1., 2., 6. };
159 
160  QGauss<elemdim> qgauss(order);
161  arma::mat map_mat = map_->element_map(ele);
162 
163  for(unsigned i=0; i<qgauss.size(); ++i) {
164  q_weights[i] = qgauss.weight(i)*weight_coefs[elemdim];
165  q_points[i] = map_->project_unit_to_real(RefElement<elemdim>::local_to_bary(qgauss.point(i)), map_mat);
166  }
167 
168  return qgauss.size();
169 }
170 
171 
172 template <int elemdim, int spacedim, class Value>
174 {
175  unsigned int ndofs = this->value_.n_rows() * this->value_.n_cols();
176  for (unsigned int k=0; k<ndofs; k++) {
177  indices[k] = (*boundary_dofs_)[ndofs*cell.idx()+k];
178  }
179  return ndofs;
180 }
181 
182 
183 template <int spacedim, class Value>
185 {
186  if (dof_indices.size() > 0)
187  WarningOut() << "Multiple initialization of FEValueHandler!";
188 
189  dh_ = init_data.dh;
190  data_vec_ = init_data.data_vec;
191  dof_indices.resize(init_data.ndofs);
192  value_.set_n_comp(init_data.n_comp);
193 }
194 
195 
196 template <int spacedim, class Value>
199 {
200  ASSERT_EQ( point_list.size(), value_list.size() ).error();
201 
202  ElementAccessor<3> cell = dh_->mesh()->element_accessor( elm.mesh_idx() ); // non-const ElementAccessor
203  if (boundary_dofs_) this->get_dof_indices( elm, dof_indices);
204  else dh_->get_dof_indices( cell, dof_indices );
205 
206  for (unsigned int k=0; k<point_list.size(); k++) {
207  Value envelope(value_list[k]);
208  envelope.zeros();
209  for (unsigned int i=0; i<dh_->fe<0>(cell)->n_dofs(); i++) {
210  envelope(i / envelope.n_cols(), i % envelope.n_rows()) += (*data_vec_)[dof_indices[i]];
211  }
212  }
213 }
214 
215 
216 template <int spacedim, class Value>
218 {
219  unsigned int ndofs = this->value_.n_rows() * this->value_.n_cols();
220  for (unsigned int k=0; k<ndofs; k++) {
221  indices[k] = (*boundary_dofs_)[ndofs*cell.idx()+k];
222  }
223  return ndofs;
224 }
225 
226 
227 template <int elemdim, int spacedim, class Value>
229 {}
230 
231 
232 // Instantiations of FEValueHandler and FEShapeHandler
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 >;
244 
245 #define INSTANCE_VALUE_HANDLER(dim) \
246 INSTANCE_VALUE_HANDLER_ALL(dim,3)
247 //INSTANCE_VALUE_HANDLER_ALL(dim,2) \
248 
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:87
RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map) const
Definition: mapping_p1.cc:291
BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map) const
Definition: mapping_p1.cc:275
VectorMPI * data_vec
Store data of Field.
void initialize(FEValueInitData init_data, MappingP1< elemdim, 3 > *map=nullptr)
Initialize data members.
VectorMPI * data_vec_
Store data of Field.
Space< spacedim >::Point Point
unsigned int mesh_idx() const
Return global idx of element in full element vector.
Definition: accessors.hh:117
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...
const FEValuesViews::Tensor< dim, spacedim > & tensor_view(unsigned int i) const
Accessor to tensor values of multicomponent FE.
Definition: fe_values.hh:398
#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:35
Symmetric Gauss-Legendre quadrature formulae on simplices.
MappingP1< elemdim, 3 > * map_
Mapping object.
double weight(const unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:161
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_
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.
Definition: fe_values.hh:388
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.
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:224
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
std::shared_ptr< std::vector< LongIdx > > boundary_dofs_
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:139
~FEValueHandler()
Destructor.
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:144
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
Calculates finite element data on the actual cell.
Definition: fe_values.hh:532
static Value::return_type fe_value(FEValues< elemdim, 3 > &fe_val, unsigned int i_dof, unsigned int i_qp)
FEValueHandler()
Constructor.
ElementMap element_map(ElementAccessor< 3 > elm) const
Definition: mapping_p1.cc:265
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.
Definition: accessors.hh:111
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
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:528