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