Flow123d  DF_patch_fe_values-dbc06cd
patch_fe_values.hh
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 patch_fe_values.hh
15  * @brief Class FEValues calculates finite element data on the actual
16  * cells such as shape function values, gradients, Jacobian of
17  * the mapping from the reference cell etc.
18  * @author Jan Stebel, David Flanderka
19  */
20 
21 #ifndef PATCH_FE_VALUES_HH_
22 #define PATCH_FE_VALUES_HH_
23 
24 
25 #include <string.h> // for memcpy
26 #include <algorithm> // for swap
27 #include <new> // for operator new[]
28 #include <string> // for operator<<
29 #include <vector> // for vector
30 #include "fem/fe_system.hh" // for FESystem
31 #include "fem/eigen_tools.hh"
33 #include "fem/patch_op.hh"
34 #include "fem/op_accessors.hh"
35 #include "mesh/ref_element.hh" // for RefElement
36 #include "mesh/accessors.hh"
38 #include "fem/arena_resource.hh"
39 #include "fem/arena_vec.hh"
40 
41 //template<unsigned int dim> class BulkValues;
42 //template<unsigned int dim> class SideValues;
43 //template<unsigned int dim> class JoinValues;
44 
45 
46 
47 template<unsigned int spacedim = 3>
49 public:
51 
53  : patch_fe_data_(1024 * 1024, 256),
56  elements_map_(300, (uint)-1)
57  {
58  element_quads_.push_back( QGauss(0, 0) );
59  for (uint dim=1; dim<4; ++dim) {
62  element_quads_.push_back( QGauss(dim, 0) );
63  }
65  }
66 
68  : PatchFEValues<spacedim>()
69  {
70  fe_ = fe;
71 
72  // TODO move initialization zero_vec_ to patch_fe_data_ constructor when we will create separate ArenaVec of DOshape functions
73  uint zero_vec_size = 300;
75  for (uint i=0; i<zero_vec_size; ++i) patch_fe_data_.zero_vec_(i) = 0.0;
76  }
77 
78 
79  /// Destructor
81  {}
82 
83  /// Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects
84  void init_finalize() {
86  }
87 
88  /// Reset PatchpointValues structures
89  void reset()
90  {
91  for (unsigned int i=0; i<spacedim; ++i) {
94  }
96  }
97 
98  /// Reinit data.
99  void reinit_patch()
100  {
101  for (auto * op : operations_) {
102  op->eval();
103  }
104  }
105 
106  /**
107  * @brief Returns the number of shape functions.
108  */
109  template<unsigned int dim>
110  inline unsigned int n_dofs() const {
111  ASSERT((dim>=0) && (dim<=3))(dim).error("Dimension must be 0, 1, 2 or 3.");
112  return fe_[Dim<dim>{}]->n_dofs();
113  }
114 
115  /**
116  * @brief Returns the number of shape functions og higher dim element.
117  */
118  template<unsigned int dim>
119  unsigned int n_dofs_high() const;
120 
121  /**
122  * @brief Returnd FiniteElement of \p component_idx for FESystem or \p fe for other types
123  */
124  template<unsigned int dim>
125  std::shared_ptr<FiniteElement<dim>> fe_comp(std::shared_ptr< FiniteElement<dim> > fe, uint component_idx) {
126  if (fe->fe_type() == FEMixedSystem) {
127  FESystem<dim> *fe_sys = dynamic_cast<FESystem<dim>*>( fe.get() );
128  return fe_sys->fe()[component_idx];
129  } else {
130  ASSERT_EQ(component_idx, 0).warning("Non-zero component_idx can only be used for FESystem.");
131  return fe;
132  }
133  }
134 
135  /// Returns pointer to FiniteElement of given dimension.
136  template<unsigned int dim>
137  std::shared_ptr<FiniteElement<dim>> fe_dim() {
138  return fe_[Dim<dim>{}];
139  }
140 
141  /** Following methods are used during update of patch. **/
142 
143  /// Clear elements_map, set values to (-1)
144  inline void clean_elements_map() {
145  std::fill(elements_map_.begin(), elements_map_.end(), (uint)-1);
146  }
147 
148  /// Add elements, sides and quadrature points registered on patch
149  template <unsigned int dim>
150  inline void add_patch_points(const DimIntegrals<dim> &integrals, ElementCacheMap *element_cache_map, std::shared_ptr<EvalPoints> eval_points) {
151  if (used_domain_[bulk_domain]) patch_point_vals_[bulk_domain][dim-1].resize_tables(eval_points->get_max_bulk_quad_size(dim), *patch_fe_data_.patch_arena_);
152  if (used_domain_[side_domain]) patch_point_vals_[side_domain][dim-1].resize_tables(eval_points->get_max_side_quad_size(dim), *patch_fe_data_.patch_arena_);
153 
154  // add bulk points
155  for (auto integral_it : integrals.bulk_) {
156  auto &patch_data = integral_it.second->patch_data();
157  for (uint i_data=0; i_data<patch_data.permanent_size(); ++i_data) {
158  uint element_patch_idx = element_cache_map->position_in_cache(patch_data[i_data].cell.elm_idx());
159  uint elm_pos = this->register_element(patch_data[i_data].cell, element_patch_idx);
160  uint i_point = 0;
161  for (auto p : integral_it.second->points(element_patch_idx) ) {
162  this->register_bulk_point(patch_data[i_data].cell, elm_pos, p.value_cache_idx(), i_point++);
163  }
164  }
165  }
166 
167  // add boundary points
168  for (auto integral_it : integrals.boundary_) {
169  auto &patch_data = integral_it.second->patch_data();
170  for (unsigned int i_data=0; i_data<patch_data.permanent_size(); ++i_data) {
171  uint side_pos = this->register_side(patch_data[i_data].side, element_cache_map);
172  uint i_point = 0;
173  for (auto p : integral_it.second->points(patch_data[i_data].side) ) {
174  this->register_side_point(patch_data[i_data].side, side_pos, p.value_cache_idx(), i_point++);
175  }
176  }
177  }
178 
179  // add edge points
180  for (auto integral_it : integrals.edge_) {
181  auto &patch_data = integral_it.second->patch_data();
182  for (unsigned int i_data=0; i_data<patch_data.permanent_size(); ++i_data) {
183  auto range = patch_data[i_data].edge_side_range;
184  for( DHCellSide edge_side : range )
185  {
186  uint side_pos = this->register_side(edge_side, element_cache_map);
187  uint i_point = 0;
188  for (auto p : integral_it.second->points(edge_side) ) {
189  this->register_side_point(edge_side, side_pos, p.value_cache_idx(), i_point++);
190  }
191  }
192  }
193  }
194 
195  // add coupling points
196  for (auto integral_it : integrals.coupling_) {
197  auto &patch_data = integral_it.second->patch_data();
198  uint side_pos, element_patch_idx, elm_pos=0;
199  uint last_element_idx = -1;
200 
201  for (unsigned int i_data=0; i_data<patch_data.permanent_size(); ++i_data) {
202  side_pos = this->register_side(patch_data[i_data].side, element_cache_map);
203  if (patch_data[i_data].cell.elm_idx() != last_element_idx) {
204  element_patch_idx = element_cache_map->position_in_cache(patch_data[i_data].cell.elm_idx());
205  elm_pos = this->register_element(patch_data[i_data].cell, element_patch_idx);
206  }
207 
208  uint i_bulk_point = 0, i_side_point = 0;
209  for (auto p_high : integral_it.second->points(patch_data[i_data].side) )
210  {
211  this->register_side_point(patch_data[i_data].side, side_pos, p_high.value_cache_idx(), i_side_point++);
212  if (patch_data[i_data].cell.elm_idx() != last_element_idx) {
213  auto p_low = p_high.lower_dim(patch_data[i_data].cell);
214  this->register_bulk_point(patch_data[i_data].cell, elm_pos, p_low.value_cache_idx(), i_bulk_point++);
215  }
216  }
217  last_element_idx = patch_data[i_data].cell.elm_idx();
218  }
219  }
220  }
221 
222  /// Register element to patch_point_vals_ table by dimension of element
223  uint register_element(DHCellAccessor cell, uint element_patch_idx) {
224  uint elem_pos = register_element_internal(cell, element_patch_idx);
226  auto map_it = ppv.n_elems_.insert( {cell.elm_idx(), ppv.i_mesh_item_} );
227  bool is_elm_added = map_it.second;
228  if (is_elm_added) {
229  ppv.int_table_(patch_elem_on_domain)(ppv.i_mesh_item_) = elem_pos;
230  ppv.i_mesh_item_++;
231  }
232  return map_it.first->second;
233  }
234 
235  /// Register side to patch_point_vals_ table by dimension of side
236  uint register_side(DHCellSide cell_side, ElementCacheMap *element_cache_map) {
237  uint dim = cell_side.dim();
238  uint element_patch_idx = element_cache_map->position_in_cache(cell_side.elem_idx());
239  uint elm_pos = register_element_internal(cell_side.cell(), element_patch_idx);
241 
242  ppv.int_table_(patch_elem_on_domain)(ppv.i_mesh_item_) = elm_pos;
243  ppv.int_table_(ref_side_on_sides)(ppv.i_mesh_item_) = cell_side.side_idx();
244  ppv.side_list_.push_back( cell_side.side() );
245  return ppv.i_mesh_item_++;
246  }
247 
248  /// Register bulk point to patch_point_vals_ table by dimension of element
249  uint register_bulk_point(DHCellAccessor cell, uint patch_elm_idx, uint elm_cache_map_idx, uint i_point_on_elem) {
250  return patch_point_vals_[bulk_domain][cell.dim()-1].register_bulk_point(patch_elm_idx, elm_cache_map_idx, cell.elm_idx(), i_point_on_elem);
251  }
252 
253  /// Register side point to patch_point_vals_ table by dimension of side
254  uint register_side_point(DHCellSide cell_side, uint patch_side_idx, uint elm_cache_map_idx, uint i_point_on_side) {
255  return patch_point_vals_[side_domain][cell_side.dim()-1].register_side_point(patch_side_idx, elm_cache_map_idx, cell_side.elem_idx(),
256  cell_side.side_idx(), i_point_on_side);
257  }
258 
259  /// return reference to assembly arena
261  return patch_fe_data_.asm_arena_;
262  }
263 
264  /// same as previous but return constant reference
265  inline const AssemblyArena &asm_arena() const {
266  return patch_fe_data_.asm_arena_;
267  }
268 
269  /// return reference to patch arena
270  inline PatchArena &patch_arena() const {
272  }
273 
274  /// Returns operation of given dim and OpType, creates it if doesn't exist
275  template<class OpType, unsigned int dim>
277  std::string op_name = typeid(OpType).name();
278  auto tpl = OperationTplHash::op_tuple(op_name, quad->size());
279  auto it = op_dependency_.find( tpl );
280  if (it == op_dependency_.end()) {
281  PatchOp<spacedim>* new_op = new OpType(*this, quad);
282  op_dependency_[tpl] = new_op;
283  operations_.push_back(new_op);
284  DebugOut().fmt("Create new operation '{}', dim: {}, quad size: {}.\n", op_name, dim, quad->size());
285  return new_op;
286  } else {
287  return it->second;
288  }
289  }
290 
291  /// Returns operation of given dim and OpType, creates it if doesn't exist
292  template<class OpType, unsigned int dim>
294  return this->template get<OpType, dim>( this->element_quad(dim) );
295  }
296 
297  /// Returns operation of given dim and OpType, creates it if doesn't exist
298  template<class OpType, unsigned int dim>
299  PatchOp<spacedim>* get(const Quadrature *quad, std::shared_ptr<FiniteElement<dim>> fe) {
300  std::string op_name = typeid(OpType).name();
301  auto tpl = OperationTplHash::op_tuple(op_name, quad->size());
302  auto it = op_dependency_.find( tpl );
303  if (it == op_dependency_.end()) {
304  PatchOp<spacedim>* new_op = new OpType(*this, quad, fe);
305  op_dependency_[tpl] = new_op;
306  operations_.push_back(new_op);
307  DebugOut().fmt("Create new operation '{}', dim: {}, quad size: {}.\n", op_name, dim, quad->size());
308  return new_op;
309  } else {
310  return it->second;
311  }
312  }
313 
314  /// Print table of all used operations - development method
315  void print_operations(ostream& stream) const {
316  stream << endl << "Table of patch FE operations:" << endl;
317  stream << std::setfill('-') << setw(160) << "" << endl;
318 
319  stream << std::setfill(' ') << " Operation" << std::setw(51) << "" << "Type" << std::setw(5) << "" << "Shape" << std::setw(2) << ""
320  << "n DOFs" << std::setw(2) << "" << "Input operations" << std::endl;
321  for (uint i=0; i<operations_.size(); ++i) {
322  stream << " " << std::left << std::setw(60) << typeid(*operations_[i]).name() << "";
323  stream << operations_[i]->dim_ << "D " << (operations_[i]->domain_ ? "side" : "bulk");
324  stream << " " << std::setw(6) << operations_[i]->format_shape() << "" << " "
325  << std::setw(7) << operations_[i]->n_dofs() << "" << " ";
326  for (auto *i_o : operations_[i]->input_ops_) stream << typeid(*i_o).name() << " ";
327  stream << std::endl;
328  }
329 
330  stream << std::setfill('=') << setw(160) << "" << endl;
331  }
332 
333  /// Getter of patch_fe_data_
335  return patch_fe_data_;
336  }
337 
338  /// Mark domain (bulk or side) as used in assembly class
339  inline void set_used_domain(fem_domain domain) {
340  used_domain_[domain] = true;
341  }
342 
343  /// Temporary method
345  ASSERT( domain<2 );
346  ASSERT( (dim>0) && (dim<=3) )(dim);
347  return patch_point_vals_[domain][dim-1];
348  }
349 
350  /// Marks data of last successfully added element to patch as permanent
352  for (uint i_dim=0; i_dim<3; ++i_dim)
353  for (uint i_domain=0; i_domain<2; ++i_domain) {
354  patch_point_vals_[i_domain][i_dim].n_mesh_items_.make_permanent();
355  }
356  }
357 
358  /// Return element quadrature (passed to element / side operations)
359  const Quadrature* element_quad(unsigned int dim) const {
360  ASSERT( dim <= 3 );
361  return &element_quads_[dim];
362  }
363 
364 private:
365  /// Register element to patch_point_vals_ table by dimension of element
368  if (elements_map_[element_patch_idx] != (uint)-1) {
369  // Return index of element on patch if it is registered repeatedly
370  return elements_map_[element_patch_idx];
371  }
372 
373  elements_map_[element_patch_idx] = ppv.elem_dim_list_->size();
374  ppv.elem_dim_list_->push_back( cell.elm() );
375  return ppv.elem_dim_list_->size() - 1;
376  }
377 
379 
380  /// Sub objects of element data of dimensions 1,2,3
382 
383  /// Sub objects of bulk and side data of dimensions 1,2,3
385 
386  MixedPtr<FiniteElement> fe_; ///< Mixed of shared pointers of FiniteElement object
387  bool used_domain_[2]; ///< Pair of flags signs holds info if bulk and side quadratures are used
388 
391 
392  /**
393  * Map of element patch indices to PatchOp::result_ and int_table_ tables
394  *
395  * TODO will be deleted after sorting elements in ElementCacheMap by dimension
396  */
398 
399  /**
400  * Array of element Quadratures of dim 0,1,2,3
401  *
402  * Items are used during construction of element/side operations. This solution solves duplicities
403  * of these operations (with different quadrature sizes). Quadrature size has no effect on result
404  * of these operations.
405  */
407 
408  friend class PatchOp<spacedim>;
409 };
410 
411 
412 
413 #endif /* PATCH_FE_VALUES_HH_ */
#define ASSERT(expr)
Definition: asserts.hh:351
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
PatchArenaResource< Resource > * get_child_arena()
Cell accessor allow iterate over DOF handler cells.
unsigned int dim() const
Return dimension of element appropriate to cell.
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Side accessor allows to iterate over sides of DOF handler cell.
unsigned int elem_idx() const
Side side() const
Return Side of given cell and side_idx.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
unsigned int dim() const
Return dimension of element appropriate to the side.
unsigned int side_idx() const
Directing class of FieldValueCache.
unsigned int position_in_cache(unsigned mesh_elm_idx, bool bdr=false) const
Return position of element stored in ElementCacheMap.
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:102
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
Definition: fe_system.hh:147
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
PatchFeData & patch_fe_data()
Getter of patch_fe_data_.
uint register_element(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
uint register_side_point(DHCellSide cell_side, uint patch_side_idx, uint elm_cache_map_idx, uint i_point_on_side)
Register side point to patch_point_vals_ table by dimension of side.
PatchPointValues< spacedim >::PatchFeData PatchFeData
std::vector< PatchOp< spacedim > * > operations_
std::vector< Quadrature > element_quads_
~PatchFEValues()
Destructor.
PatchPointValues< spacedim > & ppv(uint domain, uint dim)
Temporary method.
PatchOp< spacedim > * get_for_elem_quad()
Returns operation of given dim and OpType, creates it if doesn't exist.
PatchOp< spacedim > * get(const Quadrature *quad)
Returns operation of given dim and OpType, creates it if doesn't exist.
uint register_element_internal(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
PatchFEValues(MixedPtr< FiniteElement > fe)
void init_finalize()
Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects.
std::shared_ptr< FiniteElement< dim > > fe_dim()
Returns pointer to FiniteElement of given dimension.
OperationMap< PatchOp< spacedim > > op_dependency_
const AssemblyArena & asm_arena() const
same as previous but return constant reference
PatchArena & patch_arena() const
return reference to patch arena
uint register_bulk_point(DHCellAccessor cell, uint patch_elm_idx, uint elm_cache_map_idx, uint i_point_on_elem)
Register bulk point to patch_point_vals_ table by dimension of element.
void add_patch_points(const DimIntegrals< dim > &integrals, ElementCacheMap *element_cache_map, std::shared_ptr< EvalPoints > eval_points)
Add elements, sides and quadrature points registered on patch.
void print_operations(ostream &stream) const
Print table of all used operations - development method.
void set_used_domain(fem_domain domain)
Mark domain (bulk or side) as used in assembly class.
PatchFeData patch_fe_data_
std::vector< ElemDimList< spacedim > > elem_dim_list_vec_
Sub objects of element data of dimensions 1,2,3.
void make_permanent_ppv_data()
Marks data of last successfully added element to patch as permanent.
MixedPtr< FiniteElement > fe_
Mixed of shared pointers of FiniteElement object.
void clean_elements_map()
Clear elements_map, set values to (-1)
unsigned int n_dofs_high() const
Returns the number of shape functions og higher dim element.
std::vector< uint > elements_map_
void reinit_patch()
Reinit data.
unsigned int n_dofs() const
Returns the number of shape functions.
const Quadrature * element_quad(unsigned int dim) const
Return element quadrature (passed to element / side operations)
AssemblyArena & asm_arena()
return reference to assembly arena
bool used_domain_[2]
Pair of flags signs holds info if bulk and side quadratures are used.
uint register_side(DHCellSide cell_side, ElementCacheMap *element_cache_map)
Register side to patch_point_vals_ table by dimension of side.
PatchOp< spacedim > * get(const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Returns operation of given dim and OpType, creates it if doesn't exist.
std::vector< std::vector< PatchPointValues< spacedim > > > patch_point_vals_
Sub objects of bulk and side data of dimensions 1,2,3.
std::shared_ptr< FiniteElement< dim > > fe_comp(std::shared_ptr< FiniteElement< dim > > fe, uint component_idx)
Returnd FiniteElement of component_idx for FESystem or fe for other types.
void reset()
Reset PatchpointValues structures.
Symmetric Gauss-Legendre quadrature formulae on simplices.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
Store finite element data on the actual patch such as shape function values, gradients,...
Class FESystem for compound finite elements.
@ FEMixedSystem
std::unordered_map< std::tuple< std::string, uint >, Operation *, OperationTplHash > OperationMap
Alias for unordered_map of Operation pointer with custom hash.
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
unsigned int uint
Declares accessors to FE operations.
Base class of FE operations.
Store finite element data on the actual patch such as shape function values, gradients,...
@ patch_elem_on_domain
Index of patch element for each bulk element or side.
@ ref_side_on_sides
Ref index of side in element for each side in patch.
fem_domain
Distinguishes operations by type and size of output rows.
@ bulk_domain
@ side_domain
Definitions of particular quadrature rules on simplices.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
Set of integral of given dimension necessary in assemblation.
IntegralPtrMap< BoundaryIntegralAcc< dim > > boundary_
Boundary integrals betwwen side and boundary element of dim-1.
IntegralPtrMap< EdgeIntegralAcc< dim > > edge_
Edge integrals between elements of same dimensions.
IntegralPtrMap< BulkIntegralAcc< dim > > bulk_
Bulk integrals of elements.
IntegralPtrMap< CouplingIntegralAcc< dim > > coupling_
Coupling integrals between elements of dimensions dim and dim-1.
Definition: mixed.hh:25
static std::tuple< std::string, uint > op_tuple(std::string op_type, uint quad_size)
Create tuple from typeid(Operation).name and size of Quadrature.
AssemblyArena asm_arena_
Assembly arena, created and filled once during initialization.
ArenaVec< double > zero_vec_
ArenaVec of zero values of maximal length using in zero PatchPointValues construction.
PatchArena * patch_arena_
Patch arena, reseted before patch reinit.