Flow123d  DF_mechanic_bench-4ba6315
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  * TODO (readability, optimization):
19  * - EdgePoint::point_on is too general, can fail for non-matching sides,
20  * depends on a map;
21  * Iterator over edge sides should know index of side on the edge or
22  * we should directly iterate over pairs of sides and iterate points over both sides
23  * without mapping.
24  * - similarly for Coupling and boundary points
25  * - Points should just hold necessary indices, without reference to complex classes,
26  * these points only can be used as indices to fields to get appropriate value in the field cache
27  */
28 
29 #ifndef EVAL_SUBSET_HH_
30 #define EVAL_SUBSET_HH_
31 
32 #include <memory>
33 #include <armadillo>
34 #include "fields/eval_points.hh"
36 #include "mesh/range_wrapper.hh"
37 #include "mesh/accessors.hh"
38 #include "fem/dh_cell_accessor.hh"
39 
40 
41 class Side;
42 class BulkIntegral;
43 class EdgeIntegral;
44 class CouplingIntegral;
45 class BoundaryIntegral;
46 
47 
48 /**
49  * @brief Base point accessor class.
50  */
51 /**
52  * @brief Point accessor allow iterate over bulk quadrature points defined in local element coordinates.
53  */
54 
55 class BulkPoint {
56 public:
57  /// Default constructor
59  {}
60 
61  /// Constructor
62  BulkPoint(const ElementCacheMap *elm_cache_map, uint elem_idx, uint loc_point_idx)
63  : elm_cache_map_(elm_cache_map), elem_patch_idx_(elem_idx), local_point_idx_(loc_point_idx)
64  {}
65 
66  /// Getter of EvalPoints object.
67  inline std::shared_ptr<EvalPoints> eval_points() const {
68  ASSERT_PTR(elm_cache_map_).error("Invalid point.\n");
69  return elm_cache_map_->eval_points();
70  }
71 
72  // Getter of ElementCacheMap object.
73  inline const ElementCacheMap *elm_cache_map() const {
74  return elm_cache_map_;
75  }
76 
77 
78  // Getter of element patch index.
79  inline unsigned int elem_patch_idx() const {
80  return elem_patch_idx_;
81  }
82 
83  /// Return index in EvalPoints object
84  inline unsigned int eval_point_idx() const {
85  return local_point_idx_;
86  }
87 
88  /// Iterates to next point.
89  void inc() {
90  this->local_point_idx_++;
91  }
92 
93  /// Comparison of accessors.
94  bool operator==(const BulkPoint& other) {
96  return (local_point_idx_ == other.local_point_idx_);
97  }
98 
99 
100 protected:
102  /// Index of element in the patch.
103  unsigned int elem_patch_idx_;
104  /// Index of the local point in the integral object.
105  unsigned int local_point_idx_;
106 };
107 
108 
109 
110 
111 /**
112  * @brief General point a+ side_begin_ + ccessor allow iterate over quadrature points of given side defined in local element coordinates.
113  *
114  * Common ancestor of all side points classes (Edge-, Coupling-, BoundaryPoint)
115  */
116 class SidePoint : public BulkPoint {
117 public:
118  /// Default constructor
120  : BulkPoint() {}
121 
122  /// Constructor
124  : BulkPoint(bulk), side_begin_(side_begin)
125  {
126  //DebugOut().fmt("begin: {} sidx: {}", side_begin_, local_point_idx_);
127  }
128 
129 
130  /// Constructor
132  const EdgeIntegral *edge_integral, unsigned int local_point_idx);
133 
134  /// Return index in EvalPoints object
135  inline unsigned int eval_point_idx() const {
136  return side_begin_ + local_point_idx_;
137  }
138 
139  /// Comparison of accessors.
140  bool operator==(const SidePoint& other) {
143  return (local_point_idx_ == other.local_point_idx_);
144  }
145 
146 
147 protected:
148  //// local_point_idx_ here have meaning of index within the side.
149 
150  /// Index of side in element
151  unsigned int side_begin_;
152 
153 };
154 
155 
156 /**
157  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
158  */
159 class EdgePoint : public SidePoint {
160 public:
161  /// Default constructor
163  : SidePoint() {}
164 
165  /// Constructor
166  inline EdgePoint(BulkPoint bulk, const EdgeIntegral *edge_integral, uint side_begin);
167 
168  /// Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
169  inline EdgePoint point_on(const DHCellSide &edg_side) const;
170 
171  /// Comparison of accessors.
172  bool operator==(const EdgePoint& other) {
173  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
174  }
175 private:
177 };
178 
179 
180 /**
181  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
182  */
183 class CouplingPoint : public SidePoint {
184 public:
185  /// Default constructor
187  : SidePoint() {}
188 
189  /// Constructor
190 
191  inline CouplingPoint(BulkPoint bulk, const CouplingIntegral *coupling_integral, uint side_begin);
192 
193  /// Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
194  inline BulkPoint lower_dim(DHCellAccessor cell_lower) const;
195 
196  /// Comparison of accessors.
197  bool operator==(const CouplingPoint& other) {
198  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
199  }
200 
201 private:
202  /// Pointer to edge point set
204 };
205 
206 /**
207  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
208  */
209 class BoundaryPoint : public SidePoint {
210 public:
211  /// Default constructor
213  : SidePoint() {}
214 
215  /// Constructor
216  inline BoundaryPoint(BulkPoint bulk, const BoundaryIntegral *bdr_integral, uint side_begin);
217 
218  /// Return corresponds BulkPoint on boundary element.
219  inline BulkPoint point_bdr(ElementAccessor<3> bdr_elm) const;
220 
221  /// Comparison of accessors.
222  bool operator==(const BoundaryPoint& other) {
223  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
224  }
225 
226 private:
227  /// Pointer to edge point set
229 };
230 
231 
232 
233 /**
234  * Base integral class holds common data members and methods.
235  */
237 public:
238  /// Default constructor
239  BaseIntegral() : eval_points_(nullptr), dim_(0) {}
240 
241  /// Constructor of bulk or side subset
242  BaseIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim)
244 
245  /// Destructor
246  virtual ~BaseIntegral();
247 
248  /// Getter of eval_points
249  std::shared_ptr<EvalPoints> eval_points() const {
250  return eval_points_;
251  }
252 
253  /// Returns dimension.
254  unsigned int dim() const {
255  return dim_;
256  }
257 protected:
258  /// Pointer to EvalPoints
259  std::shared_ptr<EvalPoints> eval_points_;
260  /// Dimension of the cell on which points are placed
261  unsigned int dim_;
262 };
263 
264 /**
265  * Integral class of bulk points, allows assemblation of volume integrals.
266  */
267 class BulkIntegral : public BaseIntegral, public std::enable_shared_from_this<BulkIntegral> {
268 public:
269  /// Default constructor
271 
272  /// Constructor of bulk integral
273  BulkIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim, uint i_subset)
275  {
276  begin_idx_ = eval_points_->subset_begin(dim_, subset_index_);
277  end_idx_ = eval_points_->subset_end(dim_, subset_index_);
278  }
279 
280  /// Destructor
281  ~BulkIntegral();
282 
283  /// Return index of data block according to subset in EvalPoints object
284  inline int get_subset_idx() const {
285  return subset_index_;
286  }
287 
288 
289  /// Returns range of bulk local points for appropriate cell accessor
290  inline Range< BulkPoint > points(unsigned int element_patch_idx, const ElementCacheMap *elm_cache_map) const {
291  auto bgn_it = make_iter<BulkPoint>( BulkPoint(elm_cache_map, element_patch_idx, begin_idx_));
292  auto end_it = make_iter<BulkPoint>( BulkPoint(elm_cache_map, element_patch_idx, end_idx_));
293  return Range<BulkPoint>(bgn_it, end_it);
294  }
295 
296 protected:
297  /// Index of data block according to subset in EvalPoints object.
298  unsigned int subset_index_;
301 
302 };
303 
304 /**
305  * Integral class of side points, allows assemblation of element - element fluxes.
306  */
307 class EdgeIntegral : public BaseIntegral, public std::enable_shared_from_this<EdgeIntegral> {
308 public:
309  /// Default constructor
311  {
312  ASSERT_PERMANENT(false);
313  }
314 
315  /// Constructor of edge integral
316  EdgeIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim, uint i_subset);
317 
318  /// Destructor
319  ~EdgeIntegral();
320 
321  /// Getter of n_sides
322  inline unsigned int n_sides() const {
323  return n_sides_;
324  }
325 
326  /// Return index of data block according to subset in EvalPoints object
327  inline int get_subset_idx() const {
328  return subset_index_;
329  }
330 
331  inline uint side_begin(const DHCellSide &cell_side) const {
332  return begin_idx_ + cell_side.side_idx() * n_points_per_side_;
333  }
334 
335  /// Returns range of side local points for appropriate cell side accessor
336  inline Range< EdgePoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
337  ASSERT_EQ(cell_side.dim(), dim_);
338  //DebugOut() << "points per side: " << n_points_per_side_;
339  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
340  uint begin_idx = side_begin(cell_side);
341  auto bgn_it = make_iter<EdgePoint>( EdgePoint(
342  BulkPoint(elm_cache_map, element_patch_idx, 0), this, begin_idx));
343  auto end_it = make_iter<EdgePoint>( EdgePoint(
344  BulkPoint(elm_cache_map, element_patch_idx, n_points_per_side_), this, begin_idx));
345  return Range<EdgePoint>(bgn_it, end_it);
346  }
347 
348 
349 private:
350  unsigned int subset_index_;
352 
353  /// Number of sides (value 0 indicates bulk set)
354  unsigned int n_sides_;
355  /// Number of points. TODO: pass this to the constructor, avoid extraction from the eval_points
357 
358  friend class EvalPoints;
359  friend class EdgePoint;
360  friend class CouplingIntegral;
361  friend class BoundaryIntegral;
362 };
363 
364 /**
365  * Integral class of neighbour points, allows assemblation of element - side fluxes.
366  *
367  * Dimension corresponds with element of higher dim.
368  */
369 class CouplingIntegral : public BaseIntegral, public std::enable_shared_from_this<CouplingIntegral> {
370 public:
371  /// Default constructor
373 
374  /// Constructor of ngh integral
375  CouplingIntegral(std::shared_ptr<EdgeIntegral> edge_integral, std::shared_ptr<BulkIntegral> bulk_integral);
376 
377  /// Destructor
379 
380  /// Return index of data block according to subset of higher dim in EvalPoints object
381  inline int get_subset_high_idx() const {
382  return edge_integral_->get_subset_idx();
383  }
384 
385  /// Return index of data block according to subset of lower dim in EvalPoints object
386  inline int get_subset_low_idx() const {
387  return bulk_integral_->get_subset_idx();
388  }
389 
390  inline uint bulk_begin() const {
391  return eval_points_->subset_begin(dim_-1, bulk_integral_->get_subset_idx());
392  }
393 
394  /// Returns range of side local points for appropriate cell side accessor
395  inline Range< CouplingPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
396  ASSERT_EQ(cell_side.dim(), dim_);
397  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
398  uint begin_idx = edge_integral_->side_begin(cell_side);
399  auto bgn_it = make_iter<CouplingPoint>( CouplingPoint(
400  BulkPoint(elm_cache_map, element_patch_idx, 0), this, begin_idx) );
401  auto end_it = make_iter<CouplingPoint>( CouplingPoint(
402  BulkPoint(elm_cache_map, element_patch_idx, edge_integral_->n_points_per_side_), this, begin_idx) );;
403  return Range<CouplingPoint>(bgn_it, end_it);
404  }
405 
406 private:
407  /// Integral according to side subset part (element of higher dim) in EvalPoints object.
408  std::shared_ptr<EdgeIntegral> edge_integral_;
409  /// Integral according to bulk subset part (element of lower dim) in EvalPoints object.
410  std::shared_ptr<BulkIntegral> bulk_integral_;
411 
412  friend class CouplingPoint;
413 };
414 
415 /**
416  * Integral class of boundary points, allows assemblation of fluxes between sides and neighbouring boundary elements.
417  */
418 class BoundaryIntegral : public BaseIntegral, public std::enable_shared_from_this<BoundaryIntegral> {
419 public:
420  /// Default constructor
422 
423  /// Constructor of bulk subset
424  BoundaryIntegral(std::shared_ptr<EdgeIntegral> edge_integral, std::shared_ptr<BulkIntegral> bulk_integral);
425 
426  /// Destructor
428 
429  /// Return index of data block according to subset of higher dim in EvalPoints object
430  inline int get_subset_high_idx() const {
431  return edge_integral_->get_subset_idx();
432  }
433 
434  /// Return index of data block according to subset of lower dim (boundary) in EvalPoints object
435  inline int get_subset_low_idx() const {
436  return bulk_integral_->get_subset_idx();
437  }
438 
439  inline uint bulk_begin() const {
440  // DebugOut().fmt("edge_begin: {} bdr_begin: {}",
441  // eval_points_->subset_begin(dim_, edge_integral_->get_subset_idx()),
442  // eval_points_->subset_begin(dim_-1, bulk_integral_->get_subset_idx()));
443  return eval_points_->subset_begin(dim_-1, bulk_integral_->get_subset_idx());
444  }
445 
446  /// Returns range of bulk local points for appropriate cell accessor
447  inline Range< BoundaryPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
448  ASSERT_EQ(cell_side.dim(), dim_);
449  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
450  uint begin_idx = edge_integral_->side_begin(cell_side);
451  auto bgn_it = make_iter<BoundaryPoint>( BoundaryPoint(
452  BulkPoint(elm_cache_map, element_patch_idx, 0), this, begin_idx) );
453  auto end_it = make_iter<BoundaryPoint>( BoundaryPoint(
454  BulkPoint(elm_cache_map, element_patch_idx, edge_integral_->n_points_per_side_), this, begin_idx) );;
455  return Range<BoundaryPoint>(bgn_it, end_it);
456  }
457 
458 private:
459  /// Integral according to higher dim (bulk) element subset part in EvalPoints object.
460  std::shared_ptr<EdgeIntegral> edge_integral_;
461  /// Integral according to kower dim (boundary) element subset part in EvalPoints object.
462  std::shared_ptr<BulkIntegral> bulk_integral_;
463 
464  friend class BoundaryPoint;
465 };
466 
467 
468 /******************************************************************************
469  * Implementation of inlined methods
470  */
471 
472 /// Constructor
473 EdgePoint::EdgePoint(BulkPoint bulk, const EdgeIntegral *edge_integral, uint side_begin)
474 : SidePoint(bulk, side_begin),
475  integral_(edge_integral)
476 {}
477 
478 inline EdgePoint EdgePoint::point_on(const DHCellSide &edg_side) const {
479  uint element_patch_idx = elm_cache_map_->position_in_cache(edg_side.element().idx());
480  uint side_begin = integral_->side_begin(edg_side);
481  return EdgePoint(BulkPoint(elm_cache_map_, element_patch_idx, local_point_idx_),
482  integral_, side_begin);
483 }
484 
485 //******************************************************************************
486 CouplingPoint::CouplingPoint(BulkPoint bulk, const CouplingIntegral *coupling_integral, uint side_begin)
487 : SidePoint(bulk, side_begin),
488  integral_(coupling_integral)
489 {}
490 
491 
493  unsigned int i_elm = elm_cache_map_->position_in_cache(cell_lower.elm().idx());
494  unsigned int i_ep = integral_->bulk_begin() + local_point_idx_;
495  return BulkPoint(elm_cache_map_, i_elm, i_ep);
496 }
497 
498 
499 
500 //******************************************************************************
502 : SidePoint(bulk, side_begin),
503  integral_(bdr_integral)
504 {}
505 
506 
507 
509  unsigned int i_elm = elm_cache_map_->position_in_cache(bdr_elm.idx(), true);
510  unsigned int i_ep = integral_->bulk_begin() + local_point_idx_;
511  //DebugOut() << "begin:" << integral_->bulk_begin() << "iloc " << local_point_idx_;
512  return BulkPoint(elm_cache_map_, i_elm, i_ep);
513 }
514 
515 #endif /* EVAL_SUBSET_HH_ */
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:348
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Definition: asserts.hh:341
BaseIntegral(std::shared_ptr< EvalPoints > eval_points, unsigned int dim)
Constructor of bulk or side subset.
Definition: eval_subset.hh:242
unsigned int dim() const
Returns dimension.
Definition: eval_subset.hh:254
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points.
Definition: eval_subset.hh:249
unsigned int dim_
Dimension of the cell on which points are placed.
Definition: eval_subset.hh:261
BaseIntegral()
Default constructor.
Definition: eval_subset.hh:239
std::shared_ptr< EvalPoints > eval_points_
Pointer to EvalPoints.
Definition: eval_subset.hh:259
virtual ~BaseIntegral()
Destructor.
Definition: eval_subset.cc:28
std::shared_ptr< EdgeIntegral > edge_integral_
Integral according to higher dim (bulk) element subset part in EvalPoints object.
Definition: eval_subset.hh:460
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:447
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:435
~BoundaryIntegral()
Destructor.
Definition: eval_subset.cc:90
BoundaryIntegral()
Default constructor.
Definition: eval_subset.hh:421
std::shared_ptr< BulkIntegral > bulk_integral_
Integral according to kower dim (boundary) element subset part in EvalPoints object.
Definition: eval_subset.hh:462
friend class BoundaryPoint
Definition: eval_subset.hh:464
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
Definition: eval_subset.hh:430
uint bulk_begin() const
Definition: eval_subset.hh:439
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:209
BulkPoint point_bdr(ElementAccessor< 3 > bdr_elm) const
Return corresponds BulkPoint on boundary element.
Definition: eval_subset.hh:508
const BoundaryIntegral * integral_
Pointer to edge point set.
Definition: eval_subset.hh:228
bool operator==(const BoundaryPoint &other)
Comparison of accessors.
Definition: eval_subset.hh:222
BoundaryPoint()
Default constructor.
Definition: eval_subset.hh:212
int get_subset_idx() const
Return index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:284
BulkIntegral(std::shared_ptr< EvalPoints > eval_points, unsigned int dim, uint i_subset)
Constructor of bulk integral.
Definition: eval_subset.hh:273
~BulkIntegral()
Destructor.
Definition: eval_subset.cc:36
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:290
BulkIntegral()
Default constructor.
Definition: eval_subset.hh:270
unsigned int subset_index_
Index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:298
Base point accessor class.
Definition: eval_subset.hh:55
const ElementCacheMap * elm_cache_map_
Definition: eval_subset.hh:101
unsigned int elem_patch_idx_
Index of element in the patch.
Definition: eval_subset.hh:103
BulkPoint(const ElementCacheMap *elm_cache_map, uint elem_idx, uint loc_point_idx)
Constructor.
Definition: eval_subset.hh:62
BulkPoint()
Default constructor.
Definition: eval_subset.hh:58
std::shared_ptr< EvalPoints > eval_points() const
Getter of EvalPoints object.
Definition: eval_subset.hh:67
bool operator==(const BulkPoint &other)
Comparison of accessors.
Definition: eval_subset.hh:94
const ElementCacheMap * elm_cache_map() const
Definition: eval_subset.hh:73
unsigned int elem_patch_idx() const
Definition: eval_subset.hh:79
unsigned int local_point_idx_
Index of the local point in the integral object.
Definition: eval_subset.hh:105
void inc()
Iterates to next point.
Definition: eval_subset.hh:89
unsigned int eval_point_idx() const
Return index in EvalPoints object.
Definition: eval_subset.hh:84
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:395
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
Definition: eval_subset.hh:381
uint bulk_begin() const
Definition: eval_subset.hh:390
int get_subset_low_idx() const
Return index of data block according to subset of lower dim in EvalPoints object.
Definition: eval_subset.hh:386
CouplingIntegral()
Default constructor.
Definition: eval_subset.hh:372
std::shared_ptr< BulkIntegral > bulk_integral_
Integral according to bulk subset part (element of lower dim) in EvalPoints object.
Definition: eval_subset.hh:410
std::shared_ptr< EdgeIntegral > edge_integral_
Integral according to side subset part (element of higher dim) in EvalPoints object.
Definition: eval_subset.hh:408
~CouplingIntegral()
Destructor.
Definition: eval_subset.cc:73
friend class CouplingPoint
Definition: eval_subset.hh:412
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:183
const CouplingIntegral * integral_
Pointer to edge point set.
Definition: eval_subset.hh:203
CouplingPoint()
Default constructor.
Definition: eval_subset.hh:186
bool operator==(const CouplingPoint &other)
Comparison of accessors.
Definition: eval_subset.hh:197
BulkPoint lower_dim(DHCellAccessor cell_lower) const
Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
Definition: eval_subset.hh:492
Cell accessor allow iterate over DOF handler cells.
const 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 dim() const
Return dimension of element appropriate to the side.
ElementAccessor< 3 > element() const
unsigned int side_idx() const
int get_subset_idx() const
Return index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:327
EdgeIntegral()
Default constructor.
Definition: eval_subset.hh:310
uint n_points_per_side_
Number of points. TODO: pass this to the constructor, avoid extraction from the eval_points.
Definition: eval_subset.hh:356
unsigned int n_sides() const
Getter of n_sides.
Definition: eval_subset.hh:322
unsigned int subset_index_
Definition: eval_subset.hh:350
~EdgeIntegral()
Destructor.
Definition: eval_subset.cc:58
unsigned int n_sides_
Number of sides (value 0 indicates bulk set)
Definition: eval_subset.hh:354
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:336
friend class EdgePoint
Definition: eval_subset.hh:359
uint side_begin(const DHCellSide &cell_side) const
Definition: eval_subset.hh:331
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:159
bool operator==(const EdgePoint &other)
Comparison of accessors.
Definition: eval_subset.hh:172
EdgePoint()
Default constructor.
Definition: eval_subset.hh:162
EdgePoint point_on(const DHCellSide &edg_side) const
Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
Definition: eval_subset.hh:478
const EdgeIntegral * integral_
Definition: eval_subset.hh:176
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
Definition: accessors.hh:215
Directing class of FieldValueCache.
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points object.
unsigned int position_in_cache(unsigned mesh_elm_idx, bool bdr=false) const
Return position of element stored in ElementCacheMap.
Class holds local coordinations of evaluating points (bulk and sides) specified by element dimension.
Definition: eval_points.hh:43
Range helper class.
General point a+ side_begin_ + ccessor allow iterate over quadrature points of given side defined in ...
Definition: eval_subset.hh:116
SidePoint(BulkPoint bulk, uint side_begin)
Constructor.
Definition: eval_subset.hh:123
SidePoint()
Default constructor.
Definition: eval_subset.hh:119
unsigned int side_begin_
Index of side in element.
Definition: eval_subset.hh:151
bool operator==(const SidePoint &other)
Comparison of accessors.
Definition: eval_subset.hh:140
SidePoint(DHCellSide cell_side, const ElementCacheMap *elm_cache_map, const EdgeIntegral *edge_integral, unsigned int local_point_idx)
Constructor.
unsigned int eval_point_idx() const
Return index in EvalPoints object.
Definition: eval_subset.hh:135
@ bulk
unsigned int uint
Implementation of range helper class.