Flow123d  JS_constraints-e651b99
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 
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:
37  /**
38  * Constructor
39  *
40  * @param quad_order Order of Quadrature objects.
41  * @param asm_internals Holds shared data with GenericAssembly
42  */
43  AssemblyBasePatch(unsigned int quad_order, AssemblyInternals *asm_internals)
44  : AssemblyBasePatch<dim>() {
45  asm_internals_ = asm_internals;
46  quad_ = new QGauss(dim, 2*quad_order);
47  quad_low_ = new QGauss(dim-1, 2*quad_order);
48  }
49 
50  /// Destructor
51  virtual ~AssemblyBasePatch() {
52  delete quad_;
53  delete quad_low_;
54  }
55 
56  /**
57  * Assembles the volume integrals on cell.
58  *
59  * Method can be overridden and implemented in descendant
60  */
61  virtual inline void cell_integral(FMT_UNUSED DHCellAccessor cell, FMT_UNUSED unsigned int element_patch_idx) {}
62 
63  /**
64  * Assembles the fluxes on the boundary.
65  *
66  * Method can be overridden and implemented in descendant
67  */
68  virtual inline void boundary_side_integral(FMT_UNUSED DHCellSide cell_side) {}
69 
70  /**
71  * Assembles the fluxes between sides on the edge.
72  *
73  * Method can be overwrite and implement in descendant
74  */
75  virtual inline void edge_integral(FMT_UNUSED RangeConvert<DHEdgeSide, DHCellSide> edge_side_range) {}
76 
77  /**
78  * Assembles the fluxes between elements of different dimensions.
79  *
80  * Method can be overridden and implemented in descendant
81  */
82  virtual inline void dimjoin_intergral(FMT_UNUSED DHCellAccessor cell_lower_dim, FMT_UNUSED DHCellSide neighb_side) {}
83 
84  /**
85  * Method prepares object before assemblation (e.g. balance, ...).
86  *
87  * Method can be overridden and implemented in descendant
88  */
89  virtual void begin() {}
90 
91  /**
92  * Method finishes object after assemblation (e.g. balance, ...).
93  *
94  * Method can be overridden and implemented in descendant
95  */
96  virtual void end() {}
97 
98  /**
99  * Create and return BulkIntegral accessor of given quadrature.
100  *
101  * Method is called from descendants during construction / initialization of assembly object.
102  */
103  std::shared_ptr<BulkIntegralAcc<dim>> create_bulk_integral(Quadrature *quad) {
104  ASSERT_PERMANENT_EQ(quad->dim(), dim);
105  std::tuple<uint, uint> tpl = IntegralTplHash::integral_tuple(dim, quad->size());
106  auto result = integrals_.bulk_.insert({
107  tpl,
108  std::make_shared<BulkIntegralAcc<dim>>(asm_internals_->eval_points_, quad, &asm_internals_->fe_values_, &asm_internals_->element_cache_map_)
109  });
110  return result.first->second;
111  }
112 
113  /**
114  * Create and return EdgeIntegral accessor of given quadrature.
115  *
116  * Method is called from descendants during construction / initialization of assembly object.
117  */
118  std::shared_ptr<EdgeIntegralAcc<dim>> create_edge_integral(Quadrature *quad) {
119  ASSERT_PERMANENT_EQ(quad->dim()+1, dim);
120  std::tuple<uint, uint> tpl = IntegralTplHash::integral_tuple(dim, quad->size());
121  auto result = integrals_.edge_.insert({
122  tpl,
123  std::make_shared<EdgeIntegralAcc<dim>>(asm_internals_->eval_points_, quad, &asm_internals_->fe_values_, &asm_internals_->element_cache_map_)
124  });
125  return result.first->second;
126  }
127 
128 
129  /**
130  * Create and return CouplingIntegral accessor of given quadrature.
131  *
132  * Method is called from descendants during construction / initialization of assembly object.
133  */
134  std::shared_ptr<CouplingIntegralAcc<dim>> create_coupling_integral(Quadrature *quad) {
135  if (dim == 3) return nullptr;
136 
137  ASSERT_PERMANENT_EQ(quad->dim(), dim);
138  std::tuple<uint, uint> tpl = IntegralTplHash::integral_tuple(dim, quad->size());
139  auto result = integrals_.coupling_.insert({
140  tpl,
141  std::make_shared<CouplingIntegralAcc<dim>>(asm_internals_->eval_points_, quad, &asm_internals_->fe_values_, &asm_internals_->element_cache_map_)
142  });
143  return result.first->second;
144  }
145 
146 
147  /**
148  * Create and return BoundaryIntegral accessor of given quadrature.
149  *
150  * Method is called from descendants during construction / initialization of assembly object.
151  */
152  std::shared_ptr<BoundaryIntegralAcc<dim>> create_boundary_integral(Quadrature *quad) {
153  ASSERT_PERMANENT_EQ(quad->dim()+1, dim);
154  std::tuple<uint, uint> tpl = IntegralTplHash::integral_tuple(dim, quad->size());
155  auto result = integrals_.boundary_.insert({
156  tpl,
157  std::make_shared<BoundaryIntegralAcc<dim>>(asm_internals_->eval_points_, quad, &asm_internals_->fe_values_, &asm_internals_->element_cache_map_)
158  });
159  return result.first->second;
160  }
161 
162 
163  /**
164  * Add data of integrals to appropriate structure and register elements to ElementCacheMap.
165  *
166  * Return true if patch is full to its maximal capacity.
167  * Method is called from GenericAssembly::assembly method.
168  */
170  ASSERT_EQ(cell.dim(), dim);
171 
172  if (cell.is_own()) { // Not ghost
173  if (integrals_.bulk_.size() > 0)
174  ++( asm_internals_->fe_values_.ppv(bulk_domain, cell.dim()) ).n_mesh_items_;
175  this->add_volume_integrals(cell);
176  }
177 
178  for( DHCellSide cell_side : cell.side_range() ) {
179  if (cell.is_own()) // Not ghost
180  if ( (cell_side.side().edge().n_sides() == 1) && (cell_side.side().is_boundary()) ) {
181  this->add_boundary_integrals(cell_side);
182  }
183  if ( (cell_side.n_edge_sides() >= min_edge_sides_) && (cell_side.edge_sides().begin()->element().idx() == cell.elm_idx())) {
184  this->add_edge_integrals(cell_side);
185  }
186  }
187 
189 
191  integrals_.revert_temporary();
192  return true;
193  } else {
194  integrals_.make_permanent();
195  return false;
196  }
197  }
198 
199  /**
200  * Assembles the cell integrals for the given dimension.
201  *
202  * Method is called from GenericAssembly::assembly method.
203  */
204  virtual inline void assemble_cell_integrals() {
205  for (unsigned int i=0; i<integrals_.n_patch_cells(); ++i) {
206  this->cell_integral(integrals_.bulk_.begin()->second->patch_data()[i].cell,
207  asm_internals_->element_cache_map_.position_in_cache(integrals_.bulk_.begin()->second->patch_data()[i].cell.elm_idx())
208  );
209  }
210  // Possibly optimization but not so fast as we would assume (needs change interface of cell_integral)
211  /*for (unsigned int i=0; i<element_cache_map_->n_elements(); ++i) {
212  unsigned int elm_start = element_cache_map_->element_chunk_begin(i);
213  if (element_cache_map_->eval_point_data(elm_start).i_eval_point_ != 0) continue;
214  this->cell_integral(i, element_cache_map_->eval_point_data(elm_start).dh_loc_idx_);
215  }*/
216  }
217 
218  /**
219  * Assembles the boundary side integrals for the given dimension.
220  *
221  * Method is called from GenericAssembly::assembly method.
222  */
224  for (unsigned int i=0; i<integrals_.n_patch_boundaries(); ++i) {
225  this->boundary_side_integral( integrals_.boundary_.begin()->second->patch_data()[i].side );
226  }
227  }
228 
229  /**
230  * Assembles the edge integrals for the given dimension.
231  *
232  * Method is called from GenericAssembly::assembly method.
233  */
234  inline void assemble_edge_integrals() {
235  for (unsigned int i=0; i<integrals_.n_patch_edges(); ++i) {
236  this->edge_integral(integrals_.edge_.begin()->second->patch_data()[i].edge_side_range);
237  }
238  }
239 
240  /**
241  * Assembles the neighbours integrals for the given dimension.
242  *
243  * Method is called from GenericAssembly::assembly method.
244  */
246  for (unsigned int i=0; i<integrals_.n_patch_neighbours(); ++i) {
247  this->dimjoin_intergral(integrals_.coupling_.begin()->second->patch_data()[i].cell, integrals_.coupling_.begin()->second->patch_data()[i].side);
248  }
249  }
250 
251  /// Setter of min_edge_sides_
252  void set_min_edge_sides(unsigned int val) {
253  min_edge_sides_ = val;
254  }
255 
256  /**
257  * Clean all integral data structures
258  *
259  * Method is called from GenericAssembly::assembly method.
260  */
262  integrals_.reset();
263  }
264 
265  /// Getter of integrals_
266  const DimIntegrals<dim> &integrals() const {
267  return integrals_;
268  }
269 
270  /// Return number of DOFs
271  inline unsigned int n_dofs() {
272  return asm_internals_->fe_values_.template n_dofs<dim>();
273  }
274 
275  /// Return number of DOFs of higher dim element
276  inline unsigned int n_dofs_high() {
277  return asm_internals_->fe_values_.template n_dofs_high<dim>();
278  }
279 
280 protected:
281  /**
282  * Default constructor.
283  *
284  * Be aware if you use this constructor. Quadrature objects must be initialized manually in descendant.
285  */
287  : quad_(nullptr), quad_low_(nullptr), asm_internals_(nullptr), min_edge_sides_(2) {}
288 
289  /**
290  * Add data of volume integrals to appropriate data structure.
291  *
292  * Method is used internally in AssemblyBase
293  */
294  inline void add_volume_integrals(const DHCellAccessor &cell) {
295  for (auto integral_it : integrals_.bulk_) {
296  uint subset_idx = integral_it.second->get_subset_idx();
297  integral_it.second->patch_data().emplace_back(cell);
298 
299  unsigned int reg_idx = cell.elm().region_idx().idx();
300  // Different access than in other integrals: We can't use range method CellIntegral::points
301  // because it passes element_patch_idx as argument that is not known during patch construction.
302  for (uint i=uint( asm_internals_->eval_points_->subset_begin(dim, subset_idx) );
303  i<uint( asm_internals_->eval_points_->subset_end(dim, subset_idx) ); ++i) {
305  }
306  }
307  }
308 
309  /**
310  * Add data of edge integrals to appropriate data structure.
311  *
312  * Method is used internally in AssemblyBase
313  */
314  inline void add_edge_integrals(const DHCellSide &cell_side) {
315  auto range = cell_side.edge_sides();
316 
317  auto &ppv = asm_internals_->fe_values_.ppv(side_domain, cell_side.dim());
318  for (auto integral_it : integrals_.edge_) {
319  integral_it.second->patch_data().emplace_back(range);
320 
321  for( DHCellSide edge_side : range ) {
322  add_side_points(integral_it.second, edge_side, ppv);
323  }
324  }
325  }
326 
327  /**
328  * Add data of boundary integrals to appropriate data structure.
329  *
330  * Method is used internally in AssemblyBase
331  */
332  inline void add_boundary_integrals(const DHCellSide &bdr_side) {
333  auto &ppv = asm_internals_->fe_values_.ppv(side_domain, bdr_side.dim());
334 
335  for (auto integral_it : integrals_.boundary_) {
336  auto integral = integral_it.second;
337  integral->patch_data().emplace_back(bdr_side);
338 
339  unsigned int reg_idx = bdr_side.element().region_idx().idx();
340  ++ppv.n_mesh_items_;
341  for (auto p : integral->points(bdr_side) ) {
342  asm_internals_->element_cache_map_.add_eval_point(reg_idx, bdr_side.elem_idx(), p.eval_point_idx(), bdr_side.cell().local_idx());
343 
344  BulkPoint p_bdr = p.point_bdr(bdr_side.cond().element_accessor()); // equivalent point on boundary element
345  unsigned int bdr_reg = bdr_side.cond().element_accessor().region_idx().idx();
346  // invalid local_idx value, DHCellAccessor of boundary element doesn't exist
347  asm_internals_->element_cache_map_.add_eval_point(bdr_reg, bdr_side.cond().bc_ele_idx(), p_bdr.eval_point_idx(), -1);
348  }
349  }
350  }
351 
352  /**
353  * Add data of coupling integrals to appropriate data structure.
354  *
355  * Method is used internally in AssemblyBase
356  */
357  inline void add_coupling_integrals(const DHCellAccessor &cell) {
358  if (dim==3) return;
359  auto &ppv_low = asm_internals_->fe_values_.ppv(bulk_domain, cell.dim());
360  auto &ppv_high = asm_internals_->fe_values_.ppv(side_domain, cell.dim()+1);
361  uint i_int=0;
362  for (auto coupling_integral_it : integrals_.coupling_) {
363  auto coupling_integral = coupling_integral_it.second;
364  // Adds data of bulk points only if bulk point were not added during processing of bulk integral
365  bool add_bulk_points = !( (integrals_.bulk_.size() > 0) & cell.is_own() );
366  if (add_bulk_points) {
367  // add points of low dim element only one time and only if they have not been added in BulkIntegral
368  for( DHCellSide ngh_side : cell.neighb_sides() ) {
369  unsigned int reg_idx_low = cell.elm().region_idx().idx();
370  ++ppv_low.n_mesh_items_;
371  for (auto p : coupling_integral->points(ngh_side) ) {
372  auto p_low = p.lower_dim(cell); // equivalent point on low dim cell
373  asm_internals_->element_cache_map_.add_eval_point(reg_idx_low, cell.elm_idx(), p_low.eval_point_idx(), cell.local_idx());
374  }
375  break;
376  }
377  }
378  // Adds data of side points of all neighbour objects
379  for( DHCellSide ngh_side : cell.neighb_sides() ) { // cell -> elm lower dim, ngh_side -> elm higher dim
380  coupling_integral->patch_data().emplace_back(cell, ngh_side);
381  add_side_points(coupling_integral, ngh_side, ppv_high);
382  }
383  ++i_int;
384  }
385  }
386 
387  /**
388  * Common part of add_edge integrals and add_coupling_integrals methods
389  *
390  * Method is used internally in AssemblyBase
391  */
392  template <template <unsigned int> class IntegralAcc>
393  inline void add_side_points(std::shared_ptr< IntegralAcc<dim> > &integral, DHCellSide cell_side, PatchPointValues<3> &ppv) {
394  ++ppv.n_mesh_items_;
395  unsigned int reg_idx = cell_side.element().region_idx().idx();
396  for (auto p : integral->points(cell_side) ) {
397  asm_internals_->element_cache_map_.add_eval_point(reg_idx, cell_side.elem_idx(), p.eval_point_idx(), cell_side.cell().local_idx());
398  }
399  }
400 
401  Quadrature *quad_; ///< Quadrature used in assembling methods.
402  Quadrature *quad_low_; ///< Quadrature used in assembling methods (dim-1).
403  DimIntegrals<dim> integrals_; ///< Set of used integrals.
404  AssemblyInternals *asm_internals_; ///< Holds shared internals data with GeneriAssembly
405 
406  /**
407  * Minimal number of sides on edge.
408  *
409  * Edge integral is created and calculated if number of sides is greater or equal than this value. Default value
410  * is 2 and can be changed
411  */
412  unsigned int min_edge_sides_;
413 };
414 
415 #endif /* ASSEMBLY_BASE_HH_ */
#define ASSERT_PERMANENT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:329
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
void add_side_points(std::shared_ptr< IntegralAcc< dim > > &integral, DHCellSide cell_side, PatchPointValues< 3 > &ppv)
std::shared_ptr< CouplingIntegralAcc< dim > > create_coupling_integral(Quadrature *quad)
void assemble_neighbour_integrals()
virtual void begin()
virtual void edge_integral(FMT_UNUSED RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
unsigned int n_dofs_high()
Return number of DOFs of higher dim element.
void add_coupling_integrals(const DHCellAccessor &cell)
std::shared_ptr< BulkIntegralAcc< dim > > create_bulk_integral(Quadrature *quad)
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
bool add_integrals_of_computing_step(DHCellAccessor cell)
virtual void assemble_cell_integrals()
virtual void boundary_side_integral(FMT_UNUSED DHCellSide cell_side)
void add_boundary_integrals(const DHCellSide &bdr_side)
virtual void dimjoin_intergral(FMT_UNUSED DHCellAccessor cell_lower_dim, FMT_UNUSED DHCellSide neighb_side)
AssemblyBasePatch(unsigned int quad_order, AssemblyInternals *asm_internals)
Quadrature * quad_
Quadrature used in assembling methods.
unsigned int min_edge_sides_
AssemblyInternals * asm_internals_
Holds shared internals data with GeneriAssembly.
DimIntegrals< dim > integrals_
Set of used integrals.
virtual ~AssemblyBasePatch()
Destructor.
virtual void cell_integral(FMT_UNUSED DHCellAccessor cell, FMT_UNUSED unsigned int element_patch_idx)
std::shared_ptr< EdgeIntegralAcc< dim > > create_edge_integral(Quadrature *quad)
unsigned int n_dofs()
Return number of DOFs.
void add_edge_integrals(const DHCellSide &cell_side)
void assemble_boundary_side_integrals()
const DimIntegrals< dim > & integrals() const
Getter of integrals_.
void set_min_edge_sides(unsigned int val)
Setter of min_edge_sides_.
void add_volume_integrals(const DHCellAccessor &cell)
void assemble_edge_integrals()
virtual void end()
std::shared_ptr< BoundaryIntegralAcc< dim > > create_boundary_integral(Quadrature *quad)
uint bc_ele_idx()
Definition: accessors.hh:374
ElementAccessor< 3 > element_accessor()
Base point accessor class.
unsigned int eval_point_idx() const
Return index in EvalPoints object.
static unsigned int get()
Return number of stored elements.
Cell accessor allow iterate over DOF handler cells.
RangeConvert< DHNeighbSide, DHCellSide > neighb_sides() const
Returns range of neighbour cell of lower dimension corresponding to cell of higher dimension.
bool is_own() const
Return true if accessor represents own element (false for ghost element)
unsigned int dim() const
Return dimension of element appropriate to cell.
Range< DHCellSide > side_range() const
Returns range of cell sides.
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_.
unsigned int local_idx() const
Return local index to element (index of DOF handler).
Side accessor allows to iterate over sides of DOF handler cell.
unsigned int elem_idx() const
Boundary cond() const
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
RangeConvert< DHEdgeSide, DHCellSide > edge_sides() const
Returns range of all sides looped over common Edge.
unsigned int dim() const
Return dimension of element appropriate to the side.
ElementAccessor< 3 > element() const
RegionIdx region_idx() const
Definition: accessors.hh:201
unsigned int get_simd_rounded_size()
Returns number of eval. points with addition of max simd duplicates due to regions.
void add_eval_point(unsigned int i_reg, unsigned int i_ele, unsigned int i_eval_point, unsigned int dh_loc_idx)
unsigned int position_in_cache(unsigned mesh_elm_idx, bool bdr=false) const
Return position of element stored in ElementCacheMap.
PatchPointValues< spacedim > & ppv(uint domain, uint dim)
Temporary method.
RevertibleValue n_mesh_items_
Number of elements or sides in patch.
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
unsigned int dim() const
Definition: quadrature.hh:72
unsigned int idx() const
Returns a global index of the region.
Definition: region.hh:81
unsigned int uint
@ bulk_domain
@ side_domain
#define FMT_UNUSED
Definition: posix.h:75
Definitions of particular quadrature rules on simplices.
Holds common data shared between GenericAssemblz and Assembly<dim> classes.
PatchFEValues< 3 > fe_values_
Common FEValues object over all dimensions.
ElementCacheMap element_cache_map_
ElementCacheMap according to EvalPoints.
std::shared_ptr< EvalPoints > eval_points_
EvalPoints object shared by all integrals.
Set of integral of given dimension necessary in assemblation.
static std::tuple< uint, uint > integral_tuple(uint dim, uint quad_size)
Create tuple from dimennsion and size of Quadrature.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.