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