Flow123d  DF_patch_fe_mechanics-a6ba684
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 local index in quadrature. Temporary method - intermediate step in implementation of PatcFEValues.
135  inline unsigned int local_point_idx() const {
136  return local_point_idx_;
137  }
138 
139  /// Intermediate step in implementation of PatcFEValues.
140  virtual unsigned int side_idx() const =0;
141 
142  /// Return index in EvalPoints object
143  inline unsigned int eval_point_idx() const {
144  return side_begin_ + local_point_idx_;
145  }
146 
147  /// Comparison of accessors.
148  bool operator==(const SidePoint& other) {
151  return (local_point_idx_ == other.local_point_idx_);
152  }
153 
154 
155 protected:
156  //// local_point_idx_ here have meaning of index within the side.
157 
158  /// Index of side in element
159  unsigned int side_begin_;
160 
161 };
162 
163 
164 /**
165  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
166  */
167 class EdgePoint : public SidePoint {
168 public:
169  /// Default constructor
171  : SidePoint() {}
172 
173  /// Constructor
174  inline EdgePoint(BulkPoint bulk, const EdgeIntegral *edge_integral, uint side_begin);
175 
176  /// Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
177  inline EdgePoint point_on(const DHCellSide &edg_side) const;
178 
179  /// Intermediate step in implementation of PatcFEValues.
180  unsigned int side_idx() const override;
181 
182  /// Comparison of accessors.
183  bool operator==(const EdgePoint& other) {
184  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
185  }
186 private:
188 };
189 
190 
191 /**
192  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
193  */
194 class CouplingPoint : public SidePoint {
195 public:
196  /// Default constructor
198  : SidePoint() {}
199 
200  /// Constructor
201 
202  inline CouplingPoint(BulkPoint bulk, const CouplingIntegral *coupling_integral, uint side_begin);
203 
204  /// Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
205  inline BulkPoint lower_dim(DHCellAccessor cell_lower) const;
206 
207  /// Intermediate step in implementation of PatcFEValues.
208  unsigned int side_idx() const override;
209 
210  /// Comparison of accessors.
211  bool operator==(const CouplingPoint& other) {
212  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
213  }
214 
215 private:
216  /// Pointer to edge point set
218 };
219 
220 /**
221  * @brief Point accessor allow iterate over quadrature points of given side defined in local element coordinates.
222  */
223 class BoundaryPoint : public SidePoint {
224 public:
225  /// Default constructor
227  : SidePoint() {}
228 
229  /// Constructor
230  inline BoundaryPoint(BulkPoint bulk, const BoundaryIntegral *bdr_integral, uint side_begin);
231 
232  /// Return corresponds BulkPoint on boundary element.
233  inline BulkPoint point_bdr(ElementAccessor<3> bdr_elm) const;
234 
235  /// Intermediate step in implementation of PatcFEValues.
236  unsigned int side_idx() const override;
237 
238  /// Comparison of accessors.
239  bool operator==(const BoundaryPoint& other) {
240  return (elem_patch_idx_ == other.elem_patch_idx_) && (local_point_idx_ == other.local_point_idx_);
241  }
242 
243 private:
244  /// Pointer to edge point set
246 };
247 
248 
249 
250 /**
251  * Base integral class holds common data members and methods.
252  */
254 public:
255  /// Default constructor
256  BaseIntegral() : eval_points_(nullptr), dim_(0) {}
257 
258  /// Constructor of bulk or side subset
259  BaseIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim)
261 
262  /// Destructor
263  virtual ~BaseIntegral();
264 
265  /// Getter of eval_points
266  std::shared_ptr<EvalPoints> eval_points() const {
267  return eval_points_;
268  }
269 
270  /// Returns dimension.
271  unsigned int dim() const {
272  return dim_;
273  }
274 protected:
275  /// Pointer to EvalPoints
276  std::shared_ptr<EvalPoints> eval_points_;
277  /// Dimension of the cell on which points are placed
278  unsigned int dim_;
279 };
280 
281 /**
282  * Integral class of bulk points, allows assemblation of volume integrals.
283  */
284 class BulkIntegral : public BaseIntegral, public std::enable_shared_from_this<BulkIntegral> {
285 public:
286  /// Default constructor
288 
289  /// Constructor of bulk integral
290  BulkIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim, uint i_subset)
292  {
293  begin_idx_ = eval_points_->subset_begin(dim_, subset_index_);
294  end_idx_ = eval_points_->subset_end(dim_, subset_index_);
295  }
296 
297  /// Destructor
298  ~BulkIntegral();
299 
300  /// Return index of data block according to subset in EvalPoints object
301  inline int get_subset_idx() const {
302  return subset_index_;
303  }
304 
305 
306  /// Returns range of bulk local points for appropriate cell accessor
307  inline Range< BulkPoint > points(unsigned int element_patch_idx, const ElementCacheMap *elm_cache_map) const {
308  auto bgn_it = make_iter<BulkPoint>( BulkPoint(elm_cache_map, element_patch_idx, begin_idx_));
309  auto end_it = make_iter<BulkPoint>( BulkPoint(elm_cache_map, element_patch_idx, end_idx_));
310  return Range<BulkPoint>(bgn_it, end_it);
311  }
312 
313 protected:
314  /// Index of data block according to subset in EvalPoints object.
315  unsigned int subset_index_;
318 
319 };
320 
321 /**
322  * Integral class of side points, allows assemblation of element - element fluxes.
323  */
324 class EdgeIntegral : public BaseIntegral, public std::enable_shared_from_this<EdgeIntegral> {
325 public:
326  /// Default constructor
328  {
329  ASSERT_PERMANENT(false);
330  }
331 
332  /// Constructor of edge integral
333  EdgeIntegral(std::shared_ptr<EvalPoints> eval_points, unsigned int dim, uint i_subset);
334 
335  /// Destructor
336  ~EdgeIntegral();
337 
338  /// Getter of n_sides
339  inline unsigned int n_sides() const {
340  return n_sides_;
341  }
342 
343  /// Return index of data block according to subset in EvalPoints object
344  inline int get_subset_idx() const {
345  return subset_index_;
346  }
347 
348  inline uint side_begin(const DHCellSide &cell_side) const {
349  return begin_idx_ + cell_side.side_idx() * n_points_per_side_;
350  }
351 
352  /// Returns range of side local points for appropriate cell side accessor
353  inline Range< EdgePoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
354  ASSERT_EQ(cell_side.dim(), dim_);
355  //DebugOut() << "points per side: " << n_points_per_side_;
356  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
357  uint begin_idx = side_begin(cell_side);
358  auto bgn_it = make_iter<EdgePoint>( EdgePoint(
359  BulkPoint(elm_cache_map, element_patch_idx, 0), this, begin_idx));
360  auto end_it = make_iter<EdgePoint>( EdgePoint(
361  BulkPoint(elm_cache_map, element_patch_idx, n_points_per_side_), this, begin_idx));
362  return Range<EdgePoint>(bgn_it, end_it);
363  }
364 
365 
366 private:
367  unsigned int subset_index_;
369 
370  /// Number of sides (value 0 indicates bulk set)
371  unsigned int n_sides_;
372  /// Number of points. TODO: pass this to the constructor, avoid extraction from the eval_points
374 
375  friend class EvalPoints;
376  friend class EdgePoint;
377  friend class CouplingPoint;
378  friend class BoundaryPoint;
379  friend class CouplingIntegral;
380  friend class BoundaryIntegral;
381 };
382 
383 /**
384  * Integral class of neighbour points, allows assemblation of element - side fluxes.
385  *
386  * Dimension corresponds with element of higher dim.
387  */
388 class CouplingIntegral : public BaseIntegral, public std::enable_shared_from_this<CouplingIntegral> {
389 public:
390  /// Default constructor
392 
393  /// Constructor of ngh integral
394  CouplingIntegral(std::shared_ptr<EdgeIntegral> edge_integral, std::shared_ptr<BulkIntegral> bulk_integral);
395 
396  /// Destructor
398 
399  /// Return index of data block according to subset of higher dim in EvalPoints object
400  inline int get_subset_high_idx() const {
401  return edge_integral_->get_subset_idx();
402  }
403 
404  /// Return index of data block according to subset of lower dim in EvalPoints object
405  inline int get_subset_low_idx() const {
406  return bulk_integral_->get_subset_idx();
407  }
408 
409  inline uint bulk_begin() const {
410  return eval_points_->subset_begin(dim_-1, bulk_integral_->get_subset_idx());
411  }
412 
413  /// Returns range of side local points for appropriate cell side accessor
414  inline Range< CouplingPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
415  ASSERT_EQ(cell_side.dim(), dim_);
416  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
417  uint begin_idx = edge_integral_->side_begin(cell_side);
418  auto bgn_it = make_iter<CouplingPoint>( CouplingPoint(
419  BulkPoint(elm_cache_map, element_patch_idx, 0), this, begin_idx) );
420  auto end_it = make_iter<CouplingPoint>( CouplingPoint(
421  BulkPoint(elm_cache_map, element_patch_idx, edge_integral_->n_points_per_side_), this, begin_idx) );;
422  return Range<CouplingPoint>(bgn_it, end_it);
423  }
424 
425 private:
426  /// Integral according to side subset part (element of higher dim) in EvalPoints object.
427  std::shared_ptr<EdgeIntegral> edge_integral_;
428  /// Integral according to bulk subset part (element of lower dim) in EvalPoints object.
429  std::shared_ptr<BulkIntegral> bulk_integral_;
430 
431  friend class CouplingPoint;
432 };
433 
434 /**
435  * Integral class of boundary points, allows assemblation of fluxes between sides and neighbouring boundary elements.
436  */
437 class BoundaryIntegral : public BaseIntegral, public std::enable_shared_from_this<BoundaryIntegral> {
438 public:
439  /// Default constructor
441 
442  /// Constructor of bulk subset
443  BoundaryIntegral(std::shared_ptr<EdgeIntegral> edge_integral, std::shared_ptr<BulkIntegral> bulk_integral);
444 
445  /// Destructor
447 
448  /// Return index of data block according to subset of higher dim in EvalPoints object
449  inline int get_subset_high_idx() const {
450  return edge_integral_->get_subset_idx();
451  }
452 
453  /// Return index of data block according to subset of lower dim (boundary) in EvalPoints object
454  inline int get_subset_low_idx() const {
455  return bulk_integral_->get_subset_idx();
456  }
457 
458  inline uint bulk_begin() const {
459  // DebugOut().fmt("edge_begin: {} bdr_begin: {}",
460  // eval_points_->subset_begin(dim_, edge_integral_->get_subset_idx()),
461  // eval_points_->subset_begin(dim_-1, bulk_integral_->get_subset_idx()));
462  return eval_points_->subset_begin(dim_-1, bulk_integral_->get_subset_idx());
463  }
464 
465  /// Returns range of bulk local points for appropriate cell accessor
466  inline Range< BoundaryPoint > points(const DHCellSide &cell_side, const ElementCacheMap *elm_cache_map) const {
467  ASSERT_EQ(cell_side.dim(), dim_);
468  uint element_patch_idx = elm_cache_map->position_in_cache(cell_side.element().idx());
469  uint begin_idx = edge_integral_->side_begin(cell_side);
470  auto bgn_it = make_iter<BoundaryPoint>( BoundaryPoint(
471  BulkPoint(elm_cache_map, element_patch_idx, 0), this, begin_idx) );
472  auto end_it = make_iter<BoundaryPoint>( BoundaryPoint(
473  BulkPoint(elm_cache_map, element_patch_idx, edge_integral_->n_points_per_side_), this, begin_idx) );;
474  return Range<BoundaryPoint>(bgn_it, end_it);
475  }
476 
477 private:
478  /// Integral according to higher dim (bulk) element subset part in EvalPoints object.
479  std::shared_ptr<EdgeIntegral> edge_integral_;
480  /// Integral according to kower dim (boundary) element subset part in EvalPoints object.
481  std::shared_ptr<BulkIntegral> bulk_integral_;
482 
483  friend class BoundaryPoint;
484 };
485 
486 
487 /******************************************************************************
488  * Implementation of inlined methods
489  */
490 
491 /// Constructor
492 EdgePoint::EdgePoint(BulkPoint bulk, const EdgeIntegral *edge_integral, uint side_begin)
493 : SidePoint(bulk, side_begin),
494  integral_(edge_integral)
495 {}
496 
497 inline EdgePoint EdgePoint::point_on(const DHCellSide &edg_side) const {
498  uint element_patch_idx = elm_cache_map_->position_in_cache(edg_side.element().idx());
499  uint side_begin = integral_->side_begin(edg_side);
500  return EdgePoint(BulkPoint(elm_cache_map_, element_patch_idx, local_point_idx_),
501  integral_, side_begin);
502 }
503 
504 //******************************************************************************
505 CouplingPoint::CouplingPoint(BulkPoint bulk, const CouplingIntegral *coupling_integral, uint side_begin)
506 : SidePoint(bulk, side_begin),
507  integral_(coupling_integral)
508 {}
509 
510 
512  unsigned int i_elm = elm_cache_map_->position_in_cache(cell_lower.elm().idx());
513  unsigned int i_ep = integral_->bulk_begin() + local_point_idx_;
514  return BulkPoint(elm_cache_map_, i_elm, i_ep);
515 }
516 
517 
518 
519 //******************************************************************************
521 : SidePoint(bulk, side_begin),
522  integral_(bdr_integral)
523 {}
524 
525 
526 
528  unsigned int i_elm = elm_cache_map_->position_in_cache(bdr_elm.idx(), true);
529  unsigned int i_ep = integral_->bulk_begin() + local_point_idx_;
530  //DebugOut() << "begin:" << integral_->bulk_begin() << "iloc " << local_point_idx_;
531  return BulkPoint(elm_cache_map_, i_elm, i_ep);
532 }
533 
534 #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:259
unsigned int dim() const
Returns dimension.
Definition: eval_subset.hh:271
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points.
Definition: eval_subset.hh:266
unsigned int dim_
Dimension of the cell on which points are placed.
Definition: eval_subset.hh:278
BaseIntegral()
Default constructor.
Definition: eval_subset.hh:256
std::shared_ptr< EvalPoints > eval_points_
Pointer to EvalPoints.
Definition: eval_subset.hh:276
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:479
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:466
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:454
~BoundaryIntegral()
Destructor.
Definition: eval_subset.cc:90
BoundaryIntegral()
Default constructor.
Definition: eval_subset.hh:440
std::shared_ptr< BulkIntegral > bulk_integral_
Integral according to kower dim (boundary) element subset part in EvalPoints object.
Definition: eval_subset.hh:481
friend class BoundaryPoint
Definition: eval_subset.hh:483
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
Definition: eval_subset.hh:449
uint bulk_begin() const
Definition: eval_subset.hh:458
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:223
BulkPoint point_bdr(ElementAccessor< 3 > bdr_elm) const
Return corresponds BulkPoint on boundary element.
Definition: eval_subset.hh:527
const BoundaryIntegral * integral_
Pointer to edge point set.
Definition: eval_subset.hh:245
bool operator==(const BoundaryPoint &other)
Comparison of accessors.
Definition: eval_subset.hh:239
BoundaryPoint()
Default constructor.
Definition: eval_subset.hh:226
unsigned int side_idx() const override
Intermediate step in implementation of PatcFEValues.
Definition: eval_subset.cc:108
int get_subset_idx() const
Return index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:301
BulkIntegral(std::shared_ptr< EvalPoints > eval_points, unsigned int dim, uint i_subset)
Constructor of bulk integral.
Definition: eval_subset.hh:290
~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:307
BulkIntegral()
Default constructor.
Definition: eval_subset.hh:287
unsigned int subset_index_
Index of data block according to subset in EvalPoints object.
Definition: eval_subset.hh:315
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:414
int get_subset_high_idx() const
Return index of data block according to subset of higher dim in EvalPoints object.
Definition: eval_subset.hh:400
uint bulk_begin() const
Definition: eval_subset.hh:409
int get_subset_low_idx() const
Return index of data block according to subset of lower dim in EvalPoints object.
Definition: eval_subset.hh:405
CouplingIntegral()
Default constructor.
Definition: eval_subset.hh:391
std::shared_ptr< BulkIntegral > bulk_integral_
Integral according to bulk subset part (element of lower dim) in EvalPoints object.
Definition: eval_subset.hh:429
std::shared_ptr< EdgeIntegral > edge_integral_
Integral according to side subset part (element of higher dim) in EvalPoints object.
Definition: eval_subset.hh:427
~CouplingIntegral()
Destructor.
Definition: eval_subset.cc:73
friend class CouplingPoint
Definition: eval_subset.hh:431
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:194
unsigned int side_idx() const override
Intermediate step in implementation of PatcFEValues.
Definition: eval_subset.cc:104
const CouplingIntegral * integral_
Pointer to edge point set.
Definition: eval_subset.hh:217
CouplingPoint()
Default constructor.
Definition: eval_subset.hh:197
bool operator==(const CouplingPoint &other)
Comparison of accessors.
Definition: eval_subset.hh:211
BulkPoint lower_dim(DHCellAccessor cell_lower) const
Return corresponds EdgePoint of neighbour side of same dimension (computing of side integrals).
Definition: eval_subset.hh:511
Cell accessor allow iterate over DOF handler cells.
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:344
EdgeIntegral()
Default constructor.
Definition: eval_subset.hh:327
uint n_points_per_side_
Number of points. TODO: pass this to the constructor, avoid extraction from the eval_points.
Definition: eval_subset.hh:373
unsigned int n_sides() const
Getter of n_sides.
Definition: eval_subset.hh:339
unsigned int subset_index_
Definition: eval_subset.hh:367
~EdgeIntegral()
Destructor.
Definition: eval_subset.cc:58
unsigned int n_sides_
Number of sides (value 0 indicates bulk set)
Definition: eval_subset.hh:371
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:353
friend class EdgePoint
Definition: eval_subset.hh:376
uint side_begin(const DHCellSide &cell_side) const
Definition: eval_subset.hh:348
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
Definition: eval_subset.hh:167
bool operator==(const EdgePoint &other)
Comparison of accessors.
Definition: eval_subset.hh:183
EdgePoint()
Default constructor.
Definition: eval_subset.hh:170
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:497
unsigned int side_idx() const override
Intermediate step in implementation of PatcFEValues.
Definition: eval_subset.cc:100
const EdgeIntegral * integral_
Definition: eval_subset.hh:187
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
virtual unsigned int side_idx() const =0
Intermediate step in implementation of PatcFEValues.
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:159
bool operator==(const SidePoint &other)
Comparison of accessors.
Definition: eval_subset.hh:148
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:143
unsigned int local_point_idx() const
Return local index in quadrature. Temporary method - intermediate step in implementation of PatcFEVal...
Definition: eval_subset.hh:135
@ bulk
unsigned int uint
Implementation of range helper class.