Flow123d  JS_before_hm-1754-g1847fd3ed
eval_subset.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 eval_subset.hh
15  * @brief
16  * @author David Flanderka
17  */
18 
19 #ifndef EVAL_SUBSET_HH_
20 #define EVAL_SUBSET_HH_
21 
22 #include <memory>
23 #include <armadillo>
24 #include "fields/eval_points.hh"
26 #include "mesh/range_wrapper.hh"
27 #include "mesh/accessors.hh"
28 #include "fem/dh_cell_accessor.hh"
29 
30 
31 class Side;
32 class BulkIntegral;
33 class EdgeIntegral;
34 class CouplingIntegral;
35 class BoundaryIntegral;
36 
37 
38 /**
39  * @brief Base point accessor class.
40  *
41  * TODO: remove virtual methods from Points classes, need changes also in Field::operator().
42  */
43 class PointBase {
44 public:
45  /// Default constructor
47  : local_point_idx_(0), elm_cache_map_(nullptr) {}
48 
49  /// Constructor
50  PointBase(const ElementCacheMap *elm_cache_map, unsigned int loc_point_idx)
51  : local_point_idx_(loc_point_idx), elm_cache_map_(elm_cache_map) {}
52 
53  /// Getter of EvalPoints object.
54  inline std::shared_ptr<EvalPoints> eval_points() const {
55  ASSERT_PTR(elm_cache_map_).error("Invalid point.\n");
56  return elm_cache_map_->eval_points();
57  }
58 
59  // Getter of ElementCacheMap object.
60  inline const ElementCacheMap *elm_cache_map() const {
61  return elm_cache_map_;
62  }
63 
64  // Getter of element patch index.
65  inline unsigned int elem_patch_idx() const {
66  return elem_patch_idx_;
67  }
68 
69  /// Iterates to next point.
70  void inc() {
71  this->local_point_idx_++;
72  }
73 
74  /// Local coordinates within element
75  //template<unsigned int dim>
76  //inline arma::vec::fixed<dim> loc_coords() const {
77  // return this->eval_points()->local_point<dim>( this->eval_point_idx() );
78  //}
79 
80  // Global coordinates within element
81  //arma::vec3 coords() const;
82 
83 protected:
84  /// Index of the local point in the integral object.
85  unsigned int local_point_idx_;
86  /// Index of element in the patch.
87  unsigned int elem_patch_idx_;
88  /// Pointer ElementCacheMap needed for point evaluation.
90 };
91 
92 
93 /**
94  * @brief Point accessor allow iterate over bulk quadrature points defined in local element coordinates.
95  */
96 class BulkPoint : public PointBase {
97 public:
98  /// Default constructor
100  : PointBase() {}
101 
102  /// Constructor
104  : PointBase(elm_cache_map, cache_pos.i_ep_) {
105  this->elem_patch_idx_ = cache_pos.i_elm_;
106  }
107 
108  /// Comparison of accessors.
109  bool operator==(const BulkPoint& other) {
110  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
111  }
112 
113  /// Return index in EvalPoints object
114  inline unsigned int eval_point_idx() const {
115  return local_point_idx_;
116  }
117 };
118 
119 
120 /**
121  * @brief General point accessor allow iterate over quadrature points of given side defined in local element coordinates.
122  *
123  * Common ancestor of all side points classes (Edge-, Coupling-, BoundaryPoint)
124  */
125 class SidePoint : public PointBase {
126 public:
127  /// Default constructor
129  : PointBase() {}
130 
131  /// Constructor
132  SidePoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, unsigned int local_point_idx)
133  : PointBase(elm_cache_map, local_point_idx), side_idx_(cell_side.side_idx())
134  {
135  this->elem_patch_idx_ = this->elm_cache_map_->position_in_cache(cell_side.element().mesh_idx());
136  }
137 
138  inline unsigned int eval_point_idx() const {
139  return side_idx_ * local_point_idx_;
140  }
141 
142 protected:
143  /// Index of side in element
144  unsigned int side_idx_;
145 };
146 
147 
148 /**
149  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
150  */
151 class EdgePoint : public SidePoint {
152 public:
153  /// Default constructor
155  : SidePoint() {}
156 
157  /// Constructor
158  EdgePoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, unsigned int local_point_idx)
159  : SidePoint(cell_side, elm_cache_map, local_point_idx) {}
160 
161  /// Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
162  inline EdgePoint point_on(DHCellSide edg_side) const;
163 
164  /// Comparison of accessors.
165  bool operator==(const EdgePoint& other) {
166  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
167  }
168 };
169 
170 
171 /**
172  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
173  */
174 class CouplingPoint : public SidePoint {
175 public:
176  /// Default constructor
178  : SidePoint() {}
179 
180  /// Constructor
181  CouplingPoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, const CouplingIntegral *coupling_integral, unsigned int local_point_idx)
182  : SidePoint(cell_side, elm_cache_map, local_point_idx), integral_(coupling_integral) {}
183 
184  /// Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
185  inline BulkPoint lower_dim(DHCellAccessor cell_lower) const;
186 
187  /// Comparison of accessors.
188  bool operator==(const CouplingPoint& other) {
189  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
190  }
191 
192 private:
193  /// Pointer to edge point set
195 };
196 
197 /**
198  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
199  */
200 class BoundaryPoint : public SidePoint {
201 public:
202  /// Default constructor
204  : SidePoint() {}
205 
206  /// Constructor
207  BoundaryPoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, const BoundaryIntegral *bdr_integral, unsigned int local_point_idx)
208  : SidePoint(cell_side, elm_cache_map, local_point_idx), integral_(bdr_integral) {}
209 
210  /// Return corresponds BulkPoint on boundary element.
211  inline BulkPoint point_bdr(ElementAccessor<3> bdr_elm) const;
212 
213  /// Comparison of accessors.
214  bool operator==(const BoundaryPoint& other) {
215  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
216  }
217 
218 private:
219  /// Pointer to edge point set
221 };
222 
223 
224 
225 /**
226  * Base integral class holds common data members and methods.
227  */
229 public:
230  /// Default constructor
231  BaseIntegral() : eval_points_(nullptr), dim_(0) {}
232 
233  /// Constructor of bulk or side subset
234  BaseIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim)
236 
237  /// Destructor
238  virtual ~BaseIntegral();
239 
240  /// Getter of eval_points
241  std::shared_ptr<EvalPoints> eval_points() const {
242  return eval_points_;
243  }
244 
245  /// Returns dimension.
246  unsigned int dim() const {
247  return dim_;
248  }
249 protected:
250  /// Pointer to EvalPoints
251  std::shared_ptr<EvalPoints> eval_points_;
252  /// Dimension of points
253  unsigned int dim_;
254 };
255 
256 /**
257  * Integral class of bulk points, allows assemblation of volume integrals.
258  */
259 class BulkIntegral : public BaseIntegral, public std::enable_shared_from_this<BulkIntegral> {
260 public:
261  /// Default constructor
263 
264  /// Constructor of bulk integral
265  BulkIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim)
267 
268  /// Destructor
269  ~BulkIntegral();
270 
271  /// Return index of data block according to subset in EvalPoints object
272  inline int get_subset_idx() const {
273  return subset_index_;
274  }
275 
276  /// Returns range of bulk local points for appropriate cell accessor
277  inline Range< BulkPoint > points(unsigned int element_patch_idx, const ElementCacheMap *elm_cache_map) const {
278  auto bgn_it = make_iter<BulkPoint>( BulkPoint(elm_cache_map, PatchCacheLoc(element_patch_idx, eval_points_->subset_begin(dim_, subset_index_)) ));
279  auto end_it = make_iter<BulkPoint>( BulkPoint(elm_cache_map, PatchCacheLoc(element_patch_idx, eval_points_->subset_end(dim_, subset_index_)) ));
280  return Range<BulkPoint>(bgn_it, end_it);
281  }
282 
283 private:
284  /// Index of data block according to subset in EvalPoints object.
285  unsigned int subset_index_;
286 };
287 
288 /**
289  * Integral class of side points, allows assemblation of element - element fluxes.
290  */
291 class EdgeIntegral : public BaseIntegral, public std::enable_shared_from_this<EdgeIntegral> {
292 public:
293  /// Default constructor
295 
296  /// Constructor of edge integral
297  EdgeIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim);
298 
299  /// Destructor
300  ~EdgeIntegral();
301 
302  /// Getter of n_sides
303  inline unsigned int n_sides() const {
304  return n_sides_;
305  }
306 
307  /// Return index of data block according to subset in EvalPoints object
308  inline int get_subset_idx() const {
309  return subset_index_;
310  }
311 
312  /// Returns range of side local points for appropriate cell side accessor
313  inline Range< EdgePoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
314  unsigned int begin_idx = eval_points_->subset_begin(dim_, subset_index_);
315  unsigned int end_idx = eval_points_->subset_end(dim_, subset_index_);
316  unsigned int points_per_side = (end_idx - begin_idx) / this->n_sides();
317  auto bgn_it = make_iter<EdgePoint>( EdgePoint(cell_side, elm_cache_map, 0 ) );
318  auto end_it = make_iter<EdgePoint>( EdgePoint(cell_side, elm_cache_map, points_per_side ) );
319  return Range<EdgePoint>(bgn_it, end_it);
320  }
321 
322 
323 private:
324  /// Index of data block according to subset in EvalPoints object.
325  unsigned int subset_index_;
326  /// Number of sides (value 0 indicates bulk set)
327  unsigned int n_sides_;
328 
329  friend class EvalPoints;
330  friend class EdgePoint;
331 };
332 
333 /**
334  * Integral class of neighbour points, allows assemblation of element - side fluxes.
335  *
336  * Dimension corresponds with element of higher dim.
337  */
338 class CouplingIntegral : public BaseIntegral, public std::enable_shared_from_this<CouplingIntegral> {
339 public:
340  /// Default constructor
342 
343  /// Constructor of ngh integral
344  CouplingIntegral(std::shared_ptr<EdgeIntegral> edge_integral, std::shared_ptr<BulkIntegral> bulk_integral);
345 
346  /// Destructor
348 
349  /// Return index of data block according to subset of higher dim in EvalPoints object
350  inline int get_subset_high_idx() const {
351  return edge_integral_->get_subset_idx();
352  }
353 
354  /// Return index of data block according to subset of lower dim in EvalPoints object
355  inline int get_subset_low_idx() const {
356  return bulk_integral_->get_subset_idx();
357  }
358 
359  /// Returns range of side local points for appropriate cell side accessor
360  inline Range< CouplingPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
361  unsigned int begin_idx = eval_points_->subset_begin(dim_, get_subset_high_idx());
362  unsigned int end_idx = eval_points_->subset_end(dim_, get_subset_high_idx());
363  unsigned int points_per_side = (end_idx - begin_idx) / edge_integral_->n_sides();
364  auto bgn_it = make_iter<CouplingPoint>( CouplingPoint(cell_side, elm_cache_map, this, 0 ) );
365  auto end_it = make_iter<CouplingPoint>( CouplingPoint(cell_side, elm_cache_map, this, points_per_side ) );
366  return Range<CouplingPoint>(bgn_it, end_it);
367  }
368 
369 private:
370  /// Integral according to side subset part (element of higher dim) in EvalPoints object.
371  std::shared_ptr<EdgeIntegral> edge_integral_;
372  /// Integral according to bulk subset part (element of lower dim) in EvalPoints object.
373  std::shared_ptr<BulkIntegral> bulk_integral_;
374 
375  friend class CouplingPoint;
376 };
377 
378 /**
379  * Integral class of boundary points, allows assemblation of fluxes between sides and neighbouring boundary elements.
380  */
381 class BoundaryIntegral : public BaseIntegral, public std::enable_shared_from_this<BoundaryIntegral> {
382 public:
383  /// Default constructor
385 
386  /// Constructor of bulk subset
387  BoundaryIntegral(std::shared_ptr<EdgeIntegral> edge_integral, std::shared_ptr<BulkIntegral> bulk_integral);
388 
389  /// Destructor
391 
392  /// Return index of data block according to subset of higher dim in EvalPoints object
393  inline int get_subset_high_idx() const {
394  return edge_integral_->get_subset_idx();
395  }
396 
397  /// Return index of data block according to subset of lower dim (boundary) in EvalPoints object
398  inline int get_subset_low_idx() const {
399  return bulk_integral_->get_subset_idx();
400  }
401 
402  /// Returns range of bulk local points for appropriate cell accessor
403  inline Range< BoundaryPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
404  unsigned int begin_idx = eval_points_->subset_begin(dim_, get_subset_high_idx());
405  unsigned int end_idx = eval_points_->subset_end(dim_, get_subset_high_idx());
406  unsigned int points_per_side = (end_idx - begin_idx) / edge_integral_->n_sides();
407  auto bgn_it = make_iter<BoundaryPoint>( BoundaryPoint(cell_side, elm_cache_map, this, 0 ) );
408  auto end_it = make_iter<BoundaryPoint>( BoundaryPoint(cell_side, elm_cache_map, this, points_per_side ) );
409  return Range<BoundaryPoint>(bgn_it, end_it);
410  }
411 
412 private:
413  /// Integral according to higher dim (bulk) element subset part in EvalPoints object.
414  std::shared_ptr<EdgeIntegral> edge_integral_;
415  /// Integral according to kower dim (boundary) element subset part in EvalPoints object.
416  std::shared_ptr<BulkIntegral> bulk_integral_;
417 
418  friend class BoundaryPoint;
419 };
420 
421 
422 /******************************************************************************
423  * Implementation of inlined methods
424  */
425 
426 
427 inline EdgePoint EdgePoint::point_on(DHCellSide edg_side) const {
428  return EdgePoint(edg_side, elm_cache_map_, this->local_point_idx_);
429 }
430 
431 
433  unsigned int i_elm = elm_cache_map_->position_in_cache(cell_lower.elm().mesh_idx());
434  unsigned int i_ep = this->eval_points()->subset_begin(cell_lower.dim(), integral_->get_subset_low_idx()) + local_point_idx_;
435  PatchCacheLoc c_pos(i_elm, i_ep);
436  return BulkPoint(elm_cache_map_, c_pos);
437 }
438 
439 
440 
442  unsigned int i_elm = elm_cache_map_->position_in_cache(bdr_elm.mesh_idx());
443  unsigned int i_ep = this->eval_points()->subset_begin(bdr_elm.dim(), integral_->get_subset_low_idx()) + local_point_idx_;
444  PatchCacheLoc c_pos(i_elm, i_ep);
445  return BulkPoint(elm_cache_map_, c_pos);
446 }
447 
448 #endif /* EVAL_SUBSET_HH_ */
CouplingIntegral::points
Range< CouplingPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const
Returns range of side local points for appropriate cell side accessor.
Definition: eval_subset.hh:360
BulkIntegral::BulkIntegral
BulkIntegral()
Default constructor.
Definition: eval_subset.hh:262
ElementAccessor::dim
unsigned int dim() const
Definition: accessors.hh:157
EdgeIntegral::n_sides_
unsigned int n_sides_
Number of sides (value 0 indicates bulk set)
Definition: eval_subset.hh:327
ElementAccessor::mesh_idx
unsigned int mesh_idx() const
Return global idx of element in full element vector.
Definition: accessors.hh:187
PointBase::eval_points
std::shared_ptr< EvalPoints > eval_points() const
Getter of EvalPoints object.
Definition: eval_subset.hh:54
PointBase::elem_patch_idx_
unsigned int elem_patch_idx_
Index of element in the patch.
Definition: eval_subset.hh:87
CouplingPoint
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:174
BulkIntegral::points
Range< BulkPoint > points(unsigned int element_patch_idx, const ElementCacheMap *elm_cache_map) const
Returns range of bulk local points for appropriate cell accessor.
Definition: eval_subset.hh:277
BaseIntegral::BaseIntegral
BaseIntegral()
Default constructor.
Definition: eval_subset.hh:231
BoundaryPoint::point_bdr
BulkPoint point_bdr(ElementAccessor< 3 > bdr_elm) const
Return corresponds BulkPoint on boundary element.
Definition: eval_subset.hh:441
BoundaryIntegral::bulk_integral_
std::shared_ptr< BulkIntegral > bulk_integral_
Integral according to kower dim (boundary) element subset part in EvalPoints object.
Definition: eval_subset.hh:416
PointBase::elem_patch_idx
unsigned int elem_patch_idx() const
Definition: eval_subset.hh:65
BoundaryPoint::integral_
const BoundaryIntegral * integral_
Pointer to edge point set.
Definition: eval_subset.hh:220
CouplingIntegral::get_subset_low_idx
int get_subset_low_idx() const
Return index of data block according to subset of lower dim in EvalPoints object.
Definition: eval_subset.hh:355
CouplingIntegral::CouplingIntegral
CouplingIntegral()
Default constructor.
Definition: eval_subset.hh:341
ElementCacheMap
Directing class of FieldValueCache.
Definition: field_value_cache.hh:141
eval_points.hh
BaseIntegral::eval_points_
std::shared_ptr< EvalPoints > eval_points_
Pointer to EvalPoints.
Definition: eval_subset.hh:251
BulkPoint
Point accessor allow iterate over bulk quadrature points defined in local element coordinates.
Definition: eval_subset.hh:96
BulkIntegral::get_subset_idx
int get_subset_idx() const
Return index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:272
ElementAccessor< 3 >
PointBase
Base point accessor class.
Definition: eval_subset.hh:43
SidePoint
General point accessor allow iterate over quadrature points of given side defined in local element co...
Definition: eval_subset.hh:125
DHCellAccessor::dim
unsigned int dim() const
Return dimension of element appropriate to cell.
Definition: dh_cell_accessor.hh:101
EdgePoint
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:151
BaseIntegral::BaseIntegral
BaseIntegral(std::shared_ptr< EvalPoints > eval_points, unsigned int dim)
Constructor of bulk or side subset.
Definition: eval_subset.hh:234
ElementCacheMap::position_in_cache
unsigned int position_in_cache(unsigned mesh_elm_idx) const
Return position of element stored in ElementCacheMap.
Definition: field_value_cache.hh:205
BaseIntegral::dim_
unsigned int dim_
Dimension of points.
Definition: eval_subset.hh:253
BoundaryIntegral::get_subset_low_idx
int get_subset_low_idx() const
Return index of data block according to subset of lower dim (boundary) in EvalPoints object.
Definition: eval_subset.hh:398
BaseIntegral::dim
unsigned int dim() const
Returns dimension.
Definition: eval_subset.hh:246
PointBase::PointBase
PointBase()
Default constructor.
Definition: eval_subset.hh:46
EdgePoint::point_on
EdgePoint point_on(DHCellSide edg_side) const
Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
Definition: eval_subset.hh:427
EdgePoint::EdgePoint
EdgePoint()
Default constructor.
Definition: eval_subset.hh:154
PointBase::inc
void inc()
Iterates to next point.
Definition: eval_subset.hh:70
EdgeIntegral::~EdgeIntegral
~EdgeIntegral()
Destructor.
Definition: eval_subset.cc:49
DHCellSide
Side accessor allows to iterate over sides of DOF handler cell.
Definition: dh_cell_accessor.hh:176
BulkPoint::operator==
bool operator==(const BulkPoint &other)
Comparison of accessors.
Definition: eval_subset.hh:109
SidePoint::eval_point_idx
unsigned int eval_point_idx() const
Definition: eval_subset.hh:138
EdgeIntegral
Definition: eval_subset.hh:291
SidePoint::side_idx_
unsigned int side_idx_
Index of side in element.
Definition: eval_subset.hh:144
dh_cell_accessor.hh
accessors.hh
CouplingPoint::integral_
const CouplingIntegral * integral_
Pointer to edge point set.
Definition: eval_subset.hh:194
CouplingIntegral::CouplingPoint
friend class CouplingPoint
Definition: eval_subset.hh:375
BoundaryIntegral::BoundaryPoint
friend class BoundaryPoint
Definition: eval_subset.hh:418
BulkIntegral::subset_index_
unsigned int subset_index_
Index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:285
CouplingIntegral
Definition: eval_subset.hh:338
BoundaryIntegral::edge_integral_
std::shared_ptr< EdgeIntegral > edge_integral_
Integral according to higher dim (bulk) element subset part in EvalPoints object.
Definition: eval_subset.hh:414
EdgeIntegral::n_sides
unsigned int n_sides() const
Getter of n_sides.
Definition: eval_subset.hh:303
CouplingPoint::CouplingPoint
CouplingPoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, const CouplingIntegral *coupling_integral, unsigned int local_point_idx)
Constructor.
Definition: eval_subset.hh:181
BaseIntegral::eval_points
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points.
Definition: eval_subset.hh:241
BoundaryIntegral
Definition: eval_subset.hh:381
PatchCacheLoc
Holds pair of positions of point in cache (element and eval point)
Definition: field_value_cache.hh:81
BaseIntegral::~BaseIntegral
virtual ~BaseIntegral()
Destructor.
Definition: eval_subset.cc:28
ElementCacheMap::eval_points
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points object.
Definition: field_value_cache.hh:183
DHCellAccessor::elm
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Definition: dh_cell_accessor.hh:71
Side
Definition: accessors.hh:361
BoundaryIntegral::BoundaryIntegral
BoundaryIntegral()
Default constructor.
Definition: eval_subset.hh:384
BoundaryPoint::operator==
bool operator==(const BoundaryPoint &other)
Comparison of accessors.
Definition: eval_subset.hh:214
BulkPoint::eval_point_idx
unsigned int eval_point_idx() const
Return index in EvalPoints object.
Definition: eval_subset.hh:114
BulkIntegral::~BulkIntegral
~BulkIntegral()
Destructor.
Definition: eval_subset.cc:36
CouplingIntegral::edge_integral_
std::shared_ptr< EdgeIntegral > edge_integral_
Integral according to side subset part (element of higher dim) in EvalPoints object.
Definition: eval_subset.hh:371
EdgeIntegral::points
Range< EdgePoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const
Returns range of side local points for appropriate cell side accessor.
Definition: eval_subset.hh:313
PointBase::PointBase
PointBase(const ElementCacheMap *elm_cache_map, unsigned int loc_point_idx)
Constructor.
Definition: eval_subset.hh:50
EdgePoint::operator==
bool operator==(const EdgePoint &other)
Comparison of accessors.
Definition: eval_subset.hh:165
BoundaryIntegral::points
Range< BoundaryPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const
Returns range of bulk local points for appropriate cell accessor.
Definition: eval_subset.hh:403
DHCellSide::element
ElementAccessor< 3 > element() const
Definition: dh_cell_accessor.hh:223
PointBase::local_point_idx_
unsigned int local_point_idx_
Local coordinates within element.
Definition: eval_subset.hh:85
EdgeIntegral::get_subset_idx
int get_subset_idx() const
Return index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:308
BaseIntegral
Definition: eval_subset.hh:228
BoundaryPoint
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:200
Range
Range helper class.
Definition: range_wrapper.hh:65
CouplingPoint::lower_dim
BulkPoint lower_dim(DHCellAccessor cell_lower) const
Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
Definition: eval_subset.hh:432
BulkPoint::BulkPoint
BulkPoint(const ElementCacheMap *elm_cache_map, PatchCacheLoc cache_pos)
Constructor.
Definition: eval_subset.hh:103
DHCellAccessor
Cell accessor allow iterate over DOF handler cells.
Definition: dh_cell_accessor.hh:43
PointBase::elm_cache_map_
const ElementCacheMap * elm_cache_map_
Pointer ElementCacheMap needed for point evaluation.
Definition: eval_subset.hh:89
BoundaryPoint::BoundaryPoint
BoundaryPoint()
Default constructor.
Definition: eval_subset.hh:203
SidePoint::SidePoint
SidePoint()
Default constructor.
Definition: eval_subset.hh:128
EdgeIntegral::subset_index_
unsigned int subset_index_
Index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:325
EdgeIntegral::EdgeIntegral
EdgeIntegral()
Default constructor.
Definition: eval_subset.hh:294
PointBase::elm_cache_map
const ElementCacheMap * elm_cache_map() const
Definition: eval_subset.hh:60
field_value_cache.hh
PatchCacheLoc::i_elm_
unsigned int i_elm_
index of element in patch
Definition: field_value_cache.hh:87
CouplingPoint::CouplingPoint
CouplingPoint()
Default constructor.
Definition: eval_subset.hh:177
BoundaryIntegral::~BoundaryIntegral
~BoundaryIntegral()
Destructor.
Definition: eval_subset.cc:76
CouplingIntegral::~CouplingIntegral
~CouplingIntegral()
Destructor.
Definition: eval_subset.cc:62
SidePoint::SidePoint
SidePoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, unsigned int local_point_idx)
Constructor.
Definition: eval_subset.hh:132
BoundaryPoint::BoundaryPoint
BoundaryPoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, const BoundaryIntegral *bdr_integral, unsigned int local_point_idx)
Constructor.
Definition: eval_subset.hh:207
BoundaryIntegral::get_subset_high_idx
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
Definition: eval_subset.hh:393
BulkIntegral::BulkIntegral
BulkIntegral(std::shared_ptr< EvalPoints > eval_points, unsigned int dim)
Constructor of bulk integral.
Definition: eval_subset.hh:265
EdgePoint::EdgePoint
EdgePoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, unsigned int local_point_idx)
Constructor.
Definition: eval_subset.hh:158
BulkPoint::BulkPoint
BulkPoint()
Default constructor.
Definition: eval_subset.hh:99
ASSERT_PTR
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
CouplingPoint::operator==
bool operator==(const CouplingPoint &other)
Comparison of accessors.
Definition: eval_subset.hh:188
BulkIntegral
Definition: eval_subset.hh:259
CouplingIntegral::get_subset_high_idx
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
Definition: eval_subset.hh:350
range_wrapper.hh
Implementation of range helper class.
EdgeIntegral::EdgePoint
friend class EdgePoint
Definition: eval_subset.hh:330
CouplingIntegral::bulk_integral_
std::shared_ptr< BulkIntegral > bulk_integral_
Integral according to bulk subset part (element of lower dim) in EvalPoints object.
Definition: eval_subset.hh:373
EvalPoints
Class holds local coordinations of evaluating points (bulk and sides) specified by element dimension.
Definition: eval_points.hh:43