Flow123d  DF_patch_fe_mechanics-5faa023
assembly_base.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 assembly_base.hh
15  * @brief
16  */
17 
18 #ifndef ASSEMBLY_BASE_HH_
19 #define ASSEMBLY_BASE_HH_
20 
21 
24 #include "fem/eval_points.hh"
25 #include "fem/element_cache_map.hh"
26 #include "fem/update_flags.hh"
27 //#include "fem/op_factory.hh"
28 
29 
30 
31 /**
32  * Base class define empty methods, these methods can be overwite in descendants.
33  */
34 template <unsigned int dim>
36 {
37 public:
42 
43  /// Constructor
44  AssemblyBase(unsigned int quad_order) {
45  quad_ = new QGauss(dim, 2*quad_order);
46  quad_low_ = new QGauss(dim-1, 2*quad_order);
47  }
48 
49  /// Destructor
50  virtual ~AssemblyBase() {
51  delete quad_;
52  delete quad_low_;
53  }
54 
55  /// Assembles the volume integrals on cell.
56  virtual inline void cell_integral(FMT_UNUSED DHCellAccessor cell, FMT_UNUSED unsigned int element_patch_idx) {}
57 
58  /// Assembles the fluxes on the boundary.
59  virtual inline void boundary_side_integral(FMT_UNUSED DHCellSide cell_side) {}
60 
61  /// Assembles the fluxes between sides on the edge.
62  virtual inline void edge_integral(FMT_UNUSED RangeConvert<DHEdgeSide, DHCellSide> edge_side_range) {}
63 
64  /// Assembles the fluxes between elements of different dimensions.
65  virtual inline void dimjoin_intergral(FMT_UNUSED DHCellAccessor cell_lower_dim, FMT_UNUSED DHCellSide neighb_side) {}
66 
67  /// Method prepares object before assemblation (e.g. balance, ...).
68  virtual void begin() {}
69 
70  /// Method finishes object after assemblation (e.g. balance, ...).
71  virtual void end() {}
72 
73  /// Getter of active_integrals.
74  inline int n_active_integrals() const {
75  return active_integrals_;
76  }
77 
78  /// Create integrals according to dim of assembly object
79  void create_integrals(std::shared_ptr<EvalPoints> eval_points, AssemblyIntegrals &integrals) {
81  ASSERT_PERMANENT_PTR(quad_).error("Data member 'quad_' must be initialized if you use bulk integral!\n");
82  integrals_.bulk_ = eval_points->add_bulk<dim>(*quad_);
83  integrals.bulk_[dim-1] = integrals_.bulk_;
84  }
86  ASSERT_PERMANENT_PTR(quad_low_).error("Data member 'quad_low_' must be initialized if you use edge integral!\n");
87  integrals_.edge_ = eval_points->add_edge<dim>(*quad_low_);
88  integrals.edge_[dim-1] = integrals_.edge_;
89  }
90  if ((dim>1) && (active_integrals_ & ActiveIntegrals::coupling)) {
91  ASSERT_PERMANENT_PTR(quad_).error("Data member 'quad_' must be initialized if you use coupling integral!\n");
92  ASSERT_PERMANENT_PTR(quad_low_).error("Data member 'quad_low_' must be initialized if you use coupling integral!\n");
93  integrals_.coupling_ = eval_points->add_coupling<dim>(*quad_low_);
94  integrals.coupling_[dim-2] = integrals_.coupling_;
95  }
97  ASSERT_PERMANENT_PTR(quad_).error("Data member 'quad_' must be initialized if you use boundary integral!\n");
98  ASSERT_PERMANENT_PTR(quad_low_).error("Data member 'quad_low_' must be initialized if you use boundary integral!\n");
99  integrals_.boundary_ = eval_points->add_boundary<dim>(*quad_low_);
100  integrals.boundary_[dim-1] = integrals_.boundary_;
101  }
102  }
103 
104  /// Return BulkPoint range of appropriate dimension
105  inline Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const {
106  return integrals_.bulk_->points(element_patch_idx, element_cache_map_);
107  }
108 
109  /// Return EdgePoint range of appropriate dimension
110  inline Range< EdgePoint > edge_points(const DHCellSide &cell_side) const {
111  ASSERT( cell_side.dim() > 0 ).error("Invalid cell dimension, must be 1, 2 or 3!\n");
112  return integrals_.edge_->points(cell_side, element_cache_map_);
113  }
114 
115  /// Return CouplingPoint range of appropriate dimension
116  inline Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const {
117  ASSERT( cell_side.dim() > 1 ).error("Invalid cell dimension, must be 2 or 3!\n");
118  return integrals_.coupling_->points(cell_side, element_cache_map_);
119  }
120 
121  /// Return BoundaryPoint range of appropriate dimension
122  inline Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const {
123  ASSERT( cell_side.dim() > 0 ).error("Invalid cell dimension, must be 1, 2 or 3!\n");
124  return integrals_.boundary_->points(cell_side, element_cache_map_);
125  }
126 
127  /// Assembles the cell integrals for the given dimension.
128  virtual inline void assemble_cell_integrals(const RevertableList<BulkIntegralData> &bulk_integral_data) {
129  for (unsigned int i=0; i<bulk_integral_data.permanent_size(); ++i) {
130  if (bulk_integral_data[i].cell.dim() != dim) continue;
131  this->cell_integral(bulk_integral_data[i].cell, element_cache_map_->position_in_cache(bulk_integral_data[i].cell.elm_idx()));
132  }
133  // Possibly optimization but not so fast as we would assume (needs change interface of cell_integral)
134  /*for (unsigned int i=0; i<element_cache_map_->n_elements(); ++i) {
135  unsigned int elm_start = element_cache_map_->element_chunk_begin(i);
136  if (element_cache_map_->eval_point_data(elm_start).i_eval_point_ != 0) continue;
137  this->cell_integral(i, element_cache_map_->eval_point_data(elm_start).dh_loc_idx_);
138  }*/
139  }
140 
141  /// Assembles the boundary side integrals for the given dimension.
142  inline void assemble_boundary_side_integrals(const RevertableList<BoundaryIntegralData> &boundary_integral_data) {
143  for (unsigned int i=0; i<boundary_integral_data.permanent_size(); ++i) {
144  if (boundary_integral_data[i].side.dim() != dim) continue;
145  this->boundary_side_integral(boundary_integral_data[i].side);
146  }
147  }
148 
149  /// Assembles the edge integrals for the given dimension.
150  inline void assemble_edge_integrals(const RevertableList<EdgeIntegralData> &edge_integral_data) {
151  for (unsigned int i=0; i<edge_integral_data.permanent_size(); ++i) {
152  auto range = edge_integral_data[i].edge_side_range;
153  if (range.begin()->dim() != dim) continue;
154  this->edge_integral(edge_integral_data[i].edge_side_range);
155  }
156  }
157 
158  /// Assembles the neighbours integrals for the given dimension.
159  inline void assemble_neighbour_integrals(const RevertableList<CouplingIntegralData> &coupling_integral_data) {
160  for (unsigned int i=0; i<coupling_integral_data.permanent_size(); ++i) {
161  if (coupling_integral_data[i].side.dim() != dim) continue;
162  this->dimjoin_intergral(coupling_integral_data[i].cell, coupling_integral_data[i].side);
163  }
164  }
165 
166  /// Register cell points of volume integral
167  virtual inline void add_patch_bulk_points(FMT_UNUSED const RevertableList<BulkIntegralData> &bulk_integral_data) {}
168 
169  /// Register side points of boundary side integral
170  virtual inline void add_patch_bdr_side_points(FMT_UNUSED const RevertableList<BoundaryIntegralData> &boundary_integral_data) {}
171 
172  /// Register side points of edge integral
173  virtual inline void add_patch_edge_points(FMT_UNUSED const RevertableList<EdgeIntegralData> &edge_integral_data) {
174  }
175 
176  /// Register bulk and side points of coupling integral
177  virtual inline void add_patch_coupling_integrals(FMT_UNUSED const RevertableList<CouplingIntegralData> &coupling_integral_data) {}
178 
179 protected:
180  /// Set of integral of given dimension necessary in assemblation
181  struct DimIntegrals {
182  std::shared_ptr<BulkIntegral> bulk_; ///< Bulk integrals of elements
183  std::shared_ptr<EdgeIntegral> edge_; ///< Edge integrals between elements of same dimensions
184  std::shared_ptr<CouplingIntegral> coupling_; ///< Coupling integrals between elements of dimensions dim and dim-1
185  std::shared_ptr<BoundaryIntegral> boundary_; ///< Boundary integrals betwwen side and boundary element of dim-1
186  };
187 
188  /**
189  * Default constructor.
190  *
191  * Be aware if you use this constructor. Quadrature objects must be initialized manually in descendant.
192  */
194  : quad_(nullptr), quad_low_(nullptr) {}
195 
196  /// Print update flags to string format.
197  std::string print_update_flags(UpdateFlags u) const {
198  std::stringstream s;
199  s << u;
200  return s.str();
201  }
202 
203  Quadrature *quad_; ///< Quadrature used in assembling methods.
204  Quadrature *quad_low_; ///< Quadrature used in assembling methods (dim-1).
205  int active_integrals_; ///< Holds mask of active integrals.
206  DimIntegrals integrals_; ///< Set of used integrals.
207  ElementCacheMap *element_cache_map_; ///< ElementCacheMap shared with GenericAssembly object.
208 };
209 
210 
211 template <unsigned int dim>
212 class AssemblyBasePatch : public AssemblyBase<dim>
213 {
214 public:
219 
221  : AssemblyBase<dim>(), fe_values_(fe_values) {
222  this->quad_ = fe_values_->get_bulk_quadrature(dim);
224  }
225 
226  /// Register cell points of volume integral
227  inline void add_patch_bulk_points(const RevertableList<BulkIntegralData> &bulk_integral_data) override {
228  for (unsigned int i=0; i<bulk_integral_data.permanent_size(); ++i) {
229  if (bulk_integral_data[i].cell.dim() != dim) continue;
230  uint element_patch_idx = this->element_cache_map_->position_in_cache(bulk_integral_data[i].cell.elm_idx());
231  uint elm_pos = fe_values_->register_element(bulk_integral_data[i].cell, element_patch_idx);
232  uint i_point = 0;
233  for (auto p : this->bulk_points(element_patch_idx) ) {
234  fe_values_->register_bulk_point(bulk_integral_data[i].cell, elm_pos, p.value_cache_idx(), i_point++);
235  }
236  }
237  }
238 
239  /// Register side points of boundary side integral
240  inline void add_patch_bdr_side_points(const RevertableList<BoundaryIntegralData> &boundary_integral_data) override {
241  for (unsigned int i=0; i<boundary_integral_data.permanent_size(); ++i) {
242  if (boundary_integral_data[i].side.dim() != dim) continue;
243  uint side_pos = fe_values_->register_side(boundary_integral_data[i].side);
244  uint i_point = 0;
245  for (auto p : this->boundary_points(boundary_integral_data[i].side) ) {
246  fe_values_->register_side_point(boundary_integral_data[i].side, side_pos, p.value_cache_idx(), i_point++);
247  }
248  }
249  }
250 
251  /// Register side points of edge integral
252  inline void add_patch_edge_points(const RevertableList<EdgeIntegralData> &edge_integral_data) override {
253  for (unsigned int i=0; i<edge_integral_data.permanent_size(); ++i) {
254  auto range = edge_integral_data[i].edge_side_range;
255  if (range.begin()->dim() != dim) continue;
256  for( DHCellSide edge_side : range )
257  {
258  uint side_pos = fe_values_->register_side(edge_side);
259  uint i_point = 0;
260  for (auto p : this->edge_points(edge_side) ) {
261  fe_values_->register_side_point(edge_side, side_pos, p.value_cache_idx(), i_point++);
262  }
263  }
264  }
265  }
266 
267  /// Register bulk and side points of coupling integral
268  inline void add_patch_coupling_integrals(const RevertableList<CouplingIntegralData> &coupling_integral_data) override {
269  uint side_pos, element_patch_idx, elm_pos=0;
270  uint last_element_idx = -1;
271 
272  for (unsigned int i=0; i<coupling_integral_data.permanent_size(); ++i) {
273  if (coupling_integral_data[i].side.dim() != dim) continue;
274  side_pos = fe_values_->register_side(coupling_integral_data[i].side);
275  if (coupling_integral_data[i].cell.elm_idx() != last_element_idx) {
276  element_patch_idx = this->element_cache_map_->position_in_cache(coupling_integral_data[i].cell.elm_idx());
277  elm_pos = fe_values_->register_element(coupling_integral_data[i].cell, element_patch_idx);
278  }
279 
280  uint i_bulk_point = 0, i_side_point = 0;
281  for (auto p_high : this->coupling_points(coupling_integral_data[i].side) )
282  {
283  fe_values_->register_side_point(coupling_integral_data[i].side, side_pos, p_high.value_cache_idx(), i_side_point++);
284  if (coupling_integral_data[i].cell.elm_idx() != last_element_idx) {
285  auto p_low = p_high.lower_dim(coupling_integral_data[i].cell);
286  fe_values_->register_bulk_point(coupling_integral_data[i].cell, elm_pos, p_low.value_cache_idx(), i_bulk_point++);
287  }
288  }
289  last_element_idx = coupling_integral_data[i].cell.elm_idx();
290  }
291  }
292 
293  /// Return BulkValues object
294  inline unsigned int n_dofs() {
295  return fe_values_->template n_dofs<dim>();
296  }
297 
298  /// Return BulkValues object
300  return fe_values_->template bulk_values<dim>();
301  }
302 
303  /// Return SideValues object
305  return fe_values_->template side_values<dim>();
306  }
307 
308  /// Return JoinValues object
310  return fe_values_->template join_values<dim>();
311  }
312 
313 protected:
314  PatchFEValues<3> *fe_values_; ///< Common FEValues object over all dimensions
315 };
316 
317 
318 #endif /* ASSEMBLY_BASE_HH_ */
#define ASSERT(expr)
Definition: asserts.hh:351
#define ASSERT_PERMANENT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:337
AssemblyBasePatch(PatchFEValues< 3 > *fe_values)
void add_patch_coupling_integrals(const RevertableList< CouplingIntegralData > &coupling_integral_data) override
Register bulk and side points of coupling integral.
GenericAssemblyBase::CouplingIntegralData CouplingIntegralData
void add_patch_bulk_points(const RevertableList< BulkIntegralData > &bulk_integral_data) override
Register cell points of volume integral.
void add_patch_edge_points(const RevertableList< EdgeIntegralData > &edge_integral_data) override
Register side points of edge integral.
JoinValues< dim > join_values()
Return JoinValues object.
PatchFEValues< 3 > * fe_values_
Common FEValues object over all dimensions.
void add_patch_bdr_side_points(const RevertableList< BoundaryIntegralData > &boundary_integral_data) override
Register side points of boundary side integral.
BulkValues< dim > bulk_values()
Return BulkValues object.
GenericAssemblyBase::BoundaryIntegralData BoundaryIntegralData
unsigned int n_dofs()
Return BulkValues object.
GenericAssemblyBase::EdgeIntegralData EdgeIntegralData
GenericAssemblyBase::BulkIntegralData BulkIntegralData
SideValues< dim > side_values()
Return SideValues object.
virtual void assemble_cell_integrals(const RevertableList< BulkIntegralData > &bulk_integral_data)
Assembles the cell integrals for the given dimension.
Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const
Return BulkPoint range of appropriate dimension.
virtual void boundary_side_integral(FMT_UNUSED DHCellSide cell_side)
Assembles the fluxes on the boundary.
virtual void end()
Method finishes object after assemblation (e.g. balance, ...).
virtual void add_patch_edge_points(FMT_UNUSED const RevertableList< EdgeIntegralData > &edge_integral_data)
Register side points of edge integral.
virtual ~AssemblyBase()
Destructor.
virtual void cell_integral(FMT_UNUSED DHCellAccessor cell, FMT_UNUSED unsigned int element_patch_idx)
Assembles the volume integrals on cell.
Range< EdgePoint > edge_points(const DHCellSide &cell_side) const
Return EdgePoint range of appropriate dimension.
GenericAssemblyBase::CouplingIntegralData CouplingIntegralData
Quadrature * quad_
Quadrature used in assembling methods.
void create_integrals(std::shared_ptr< EvalPoints > eval_points, AssemblyIntegrals &integrals)
Create integrals according to dim of assembly object.
void assemble_edge_integrals(const RevertableList< EdgeIntegralData > &edge_integral_data)
Assembles the edge integrals for the given dimension.
virtual void add_patch_coupling_integrals(FMT_UNUSED const RevertableList< CouplingIntegralData > &coupling_integral_data)
Register bulk and side points of coupling integral.
DimIntegrals integrals_
Set of used integrals.
virtual void begin()
Method prepares object before assemblation (e.g. balance, ...).
std::string print_update_flags(UpdateFlags u) const
Print update flags to string format.
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
virtual void add_patch_bdr_side_points(FMT_UNUSED const RevertableList< BoundaryIntegralData > &boundary_integral_data)
Register side points of boundary side integral.
GenericAssemblyBase::BulkIntegralData BulkIntegralData
void assemble_neighbour_integrals(const RevertableList< CouplingIntegralData > &coupling_integral_data)
Assembles the neighbours integrals for the given dimension.
Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const
Return CouplingPoint range of appropriate dimension.
virtual void dimjoin_intergral(FMT_UNUSED DHCellAccessor cell_lower_dim, FMT_UNUSED DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
virtual void edge_integral(FMT_UNUSED RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Assembles the fluxes between sides on the edge.
int active_integrals_
Holds mask of active integrals.
Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const
Return BoundaryPoint range of appropriate dimension.
AssemblyBase(unsigned int quad_order)
Constructor.
void assemble_boundary_side_integrals(const RevertableList< BoundaryIntegralData > &boundary_integral_data)
Assembles the boundary side integrals for the given dimension.
GenericAssemblyBase::EdgeIntegralData EdgeIntegralData
virtual void add_patch_bulk_points(FMT_UNUSED const RevertableList< BulkIntegralData > &bulk_integral_data)
Register cell points of volume integral.
int n_active_integrals() const
Getter of active_integrals.
ElementCacheMap * element_cache_map_
ElementCacheMap shared with GenericAssembly object.
GenericAssemblyBase::BoundaryIntegralData BoundaryIntegralData
Cell accessor allow iterate over DOF handler cells.
Side accessor allows to iterate over sides of DOF handler cell.
unsigned int dim() const
Return dimension of element appropriate to the side.
Directing class of FieldValueCache.
unsigned int position_in_cache(unsigned mesh_elm_idx, bool bdr=false) const
Return position of element stored in ElementCacheMap.
uint register_element(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx, uint i_point_on_elem)
Register bulk point to patch_point_vals_ table by dimension of element.
Quadrature * get_bulk_quadrature(uint dim) const
Getter for bulk quadrature of given dimension.
Quadrature * get_side_quadrature(uint dim) const
Getter for side quadrature of given dimension.
uint register_side(DHCellSide cell_side)
Register side to patch_point_vals_ table by dimension of side.
uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx, uint i_point_on_side)
Register side point to patch_point_vals_ table by dimension of side.
Symmetric Gauss-Legendre quadrature formulae on simplices.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Range helper class.
@ coupling
@ boundary
@ bulk
@ edge
unsigned int uint
#define FMT_UNUSED
Definition: posix.h:75
Definitions of particular quadrature rules on simplices.
Set of integral of given dimension necessary in assemblation.
std::shared_ptr< EdgeIntegral > edge_
Edge integrals between elements of same dimensions.
std::shared_ptr< CouplingIntegral > coupling_
Coupling integrals between elements of dimensions dim and dim-1.
std::shared_ptr< BulkIntegral > bulk_
Bulk integrals of elements.
std::shared_ptr< BoundaryIntegral > boundary_
Boundary integrals betwwen side and boundary element of dim-1.
Set of all used integral necessary in assemblation.
std::array< std::shared_ptr< BulkIntegral >, 3 > bulk_
Bulk integrals of elements of dimensions 1, 2, 3.
std::array< std::shared_ptr< CouplingIntegral >, 2 > coupling_
Coupling integrals between elements of dimensions 1-2, 2-3.
std::array< std::shared_ptr< BoundaryIntegral >, 3 > boundary_
Boundary integrals betwwen elements of dimensions 1, 2, 3 and boundaries.
std::array< std::shared_ptr< EdgeIntegral >, 3 > edge_
Edge integrals between elements of dimensions 1, 2, 3.
std::size_t permanent_size() const
Return permanent size of list.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68