Flow123d  DF_patch_fe_darcy_complete-579fe1e
eval_points.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_points.hh
15  * @brief
16  * @author David Flanderka
17  */
18 
19 #ifndef EVAL_POINTS_HH_
20 #define EVAL_POINTS_HH_
21 
22 
23 #include <vector>
24 #include <memory>
25 #include <armadillo>
26 #include "fem/integral_data.hh"
27 #include "mesh/range_wrapper.hh"
28 #include "system/asserts.hh"
29 #include "system/armor.hh"
30 
31 class Side;
32 class Quadrature;
33 class ElementCacheMap;
34 class BulkIntegral;
35 class EdgeIntegral;
36 class CouplingIntegral;
37 class BoundaryIntegral;
38 template <unsigned int dim> class BulkIntegralAcc;
39 template <unsigned int dim> class EdgeIntegralAcc;
40 template <unsigned int dim> class CouplingIntegralAcc;
41 template <unsigned int dim> class BoundaryIntegralAcc;
42 template <int spacedim> class ElementAccessor;
43 template <unsigned int spacedim> class PatchFEValues;
44 namespace internal_integrals {
45  class Bulk;
46  class Edge;
47 }
48 
49 
50 /**
51  * @brief Class holds local coordinations of evaluating points (bulk and sides)
52  * specified by element dimension.
53  */
54 class EvalPoints : public std::enable_shared_from_this<EvalPoints> {
55 public:
56  /// Undefined dimension of new (empty) object
57  static const unsigned int undefined_dim;
58 
59  /// Maximal number of hold subsets.
60  static constexpr unsigned int max_subsets = 10;
61 
62  /// Maximal average number of points hold in subset.
63  static const unsigned int max_subset_points = 30;
64 
65  /// Constructor
66  EvalPoints();
67 
68  /// Return size of evaluation points object (number of points).
69  inline unsigned int size(unsigned int dim) const {
70  return dim_eval_points_[dim].size();
71  }
72 
73  /// Return local coordinates of given local point and appropriate dim.
74  template<unsigned int dim>
75  inline arma::vec::fixed<dim> local_point(unsigned int local_point_idx) const {
76  ASSERT_GT(dim, 0).error("Dimension 0 not supported!\n");
77  return dim_eval_points_[dim].local_point<dim>(local_point_idx);
78  }
79 
80  /// Return begin index of appropriate subset data.
81  inline int subset_begin(unsigned int dim, unsigned int idx) const {
82  return dim_eval_points_[dim].subset_begin(idx);
83  }
84 
85  /// Return end index of appropriate subset data.
86  inline int subset_end(unsigned int dim, unsigned int idx) const {
87  return dim_eval_points_[dim].subset_end(idx);
88  }
89 
90  /// Return number of local points corresponding to subset.
91  inline int subset_size(unsigned int dim, unsigned int idx) const {
92  return dim_eval_points_[dim].subset_size(idx);
93  }
94 
95  /// Return number of subsets.
96  inline unsigned int n_subsets(unsigned int dim) const {
97  return dim_eval_points_[dim].n_subsets();
98  }
99 
100  /**
101  * Registers point set from quadrature.
102  *
103  * Returns an object referencing to the EvalPoints and list of its points.
104  */
105  template <unsigned int dim>
106  std::shared_ptr<internal_integrals::Bulk> add_bulk_internal(Quadrature *);
107 
108  /// The same as add_bulk but for edge points on sides.
109  template <unsigned int dim>
110  std::shared_ptr<internal_integrals::Edge> add_edge_internal(Quadrature *);
111 
112  /// Return maximal size of evaluation points objects.
113  inline unsigned int max_size() const {
114  return max_size_;
115  }
116 
117  inline void clear() {
118  for (uint i=0; i<4; ++i)
119  dim_eval_points_[i].clear();
120  }
121 
122  /// Return maximal size of quadrature of given dimension of bulk integral (Bulk, Coupling /lower-dim/)
123  uint get_max_bulk_quad_size(unsigned int dim) const;
124 
125  /// Return maximal size of quadrature of given dimension of side integral (Edge, Coupling /higher-dim/, Boundary)
126  uint get_max_side_quad_size(unsigned int dim) const;
127 
128 private:
129  /// Subobject holds evaluation points data of one dimension (0,1,2,3)
131  public:
132  /// Constructor
133  DimEvalPoints(unsigned int dim);
134 
135  /// Return size of evaluation points object (number of points).
136  inline unsigned int size() const {
137  if (dim_==0) return n_subsets_;
138  else return local_points_.size();
139  }
140 
141  /// Return local coordinates of given local point.
142  template<unsigned int dim>
143  inline arma::vec::fixed<dim> local_point(unsigned int local_point_idx) const {
144  ASSERT_LT(local_point_idx, this->size());
145  return local_points_.vec<dim>(local_point_idx);
146  }
147 
148  /// Return begin index of appropriate subset data.
149  inline int subset_begin(unsigned int idx) const {
150  ASSERT_LT(idx, n_subsets());
151  return subset_starts_[idx];
152  }
153 
154  /// Return end index of appropriate subset data.
155  inline int subset_end(unsigned int idx) const {
156  ASSERT_LT(idx, n_subsets());
157  return subset_starts_[idx+1];
158  }
159 
160  /// Return number of local points corresponding to subset.
161  inline int subset_size(unsigned int idx) const {
162  ASSERT_LT(idx, n_subsets());
163  return subset_starts_[idx+1] - subset_starts_[idx];
164  }
165 
166  /// Return number of subsets.
167  inline unsigned int n_subsets() const {
168  return n_subsets_;
169  }
170 
171  /// Adds set of local point to local_points_ (bulk or side).
172  template <unsigned int dim>
173  void add_local_points(const Armor::Array<double> & quad_points);
174 
175  /// Adds new subset and its end size to subset_starts_ array.
176  uint add_subset();
177 
178  inline void clear() {
180  n_subsets_ = 0;
181  }
182  private:
183  Armor::Array<double> local_points_; ///< Local coords of points vector
184  std::array<int, EvalPoints::max_subsets+1> subset_starts_; ///< Indices of subsets data in local_points_ vector, used size is n_subsets_ + 1
185  unsigned int n_subsets_; ///< Number of subset
186  unsigned int dim_; ///< Dimension of local points
187 
188  friend class EvalPoints;
189  };
190 
191  inline void set_max_size() {
192  max_size_ = std::max( std::max( size(0), size(1) ), std::max( size(2), size(3) ) );
193  }
194 
195  /// Common implementation of get_max_bulk_quad_size and get_max_side_quad_size
196  template<class Integral>
197  uint get_max_integral_quad_size(IntegralPtrMap<Integral> integrals, unsigned int dim) const;
198 
199 
200  /// Sub objects of dimensions 0,1,2,3
201  std::array<DimEvalPoints, 4> dim_eval_points_;
202 
203  /// Maps of all BulkIntegrals of dimensions 0,1,2,3
205 
206  /// Maps of all EdgeIntegrals of dimensions 1,2,3
208 
209  /// Maximal number of used EvalPoints.
210  unsigned int max_size_;
211 
212 };
213 
214 
215 #endif /* EVAL_POINTS_HH_ */
Definitions of ASSERTS.
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:301
#define ASSERT_GT(a, b)
Definition of comparative assert macro (Greater Than) only for debug mode.
Definition: asserts.hh:317
void resize(uint size)
Definition: armor.hh:758
ArmaVec< Type, nr > vec(uint mat_index) const
Definition: armor.hh:869
unsigned int size() const
Definition: armor.hh:776
Directing class of FieldValueCache.
Subobject holds evaluation points data of one dimension (0,1,2,3)
Definition: eval_points.hh:130
unsigned int n_subsets_
Number of subset.
Definition: eval_points.hh:185
int subset_begin(unsigned int idx) const
Return begin index of appropriate subset data.
Definition: eval_points.hh:149
std::array< int, EvalPoints::max_subsets+1 > subset_starts_
Indices of subsets data in local_points_ vector, used size is n_subsets_ + 1.
Definition: eval_points.hh:184
unsigned int size() const
Return size of evaluation points object (number of points).
Definition: eval_points.hh:136
arma::vec::fixed< dim > local_point(unsigned int local_point_idx) const
Return local coordinates of given local point.
Definition: eval_points.hh:143
DimEvalPoints(unsigned int dim)
Constructor.
Definition: eval_points.cc:107
unsigned int n_subsets() const
Return number of subsets.
Definition: eval_points.hh:167
uint add_subset()
Adds new subset and its end size to subset_starts_ array.
Definition: eval_points.cc:127
unsigned int dim_
Dimension of local points.
Definition: eval_points.hh:186
int subset_end(unsigned int idx) const
Return end index of appropriate subset data.
Definition: eval_points.hh:155
Armor::Array< double > local_points_
Local coords of points vector.
Definition: eval_points.hh:183
void add_local_points(const Armor::Array< double > &quad_points)
Adds set of local point to local_points_ (bulk or side).
Definition: eval_points.cc:115
int subset_size(unsigned int idx) const
Return number of local points corresponding to subset.
Definition: eval_points.hh:161
Class holds local coordinations of evaluating points (bulk and sides) specified by element dimension.
Definition: eval_points.hh:54
int subset_begin(unsigned int dim, unsigned int idx) const
Return begin index of appropriate subset data.
Definition: eval_points.hh:81
void set_max_size()
Definition: eval_points.hh:191
static const unsigned int undefined_dim
Undefined dimension of new (empty) object.
Definition: eval_points.hh:57
unsigned int n_subsets(unsigned int dim) const
Return number of subsets.
Definition: eval_points.hh:96
unsigned int size(unsigned int dim) const
Return size of evaluation points object (number of points).
Definition: eval_points.hh:69
std::shared_ptr< internal_integrals::Edge > add_edge_internal(Quadrature *)
The same as add_bulk but for edge points on sides.
Definition: eval_points.cc:69
unsigned int max_size() const
Return maximal size of evaluation points objects.
Definition: eval_points.hh:113
int subset_end(unsigned int dim, unsigned int idx) const
Return end index of appropriate subset data.
Definition: eval_points.hh:86
std::shared_ptr< internal_integrals::Bulk > add_bulk_internal(Quadrature *)
Definition: eval_points.cc:35
unsigned int max_size_
Maximal number of used EvalPoints.
Definition: eval_points.hh:210
static const unsigned int max_subset_points
Maximal average number of points hold in subset.
Definition: eval_points.hh:63
uint get_max_integral_quad_size(IntegralPtrMap< Integral > integrals, unsigned int dim) const
Common implementation of get_max_bulk_quad_size and get_max_side_quad_size.
Definition: eval_points.cc:98
void clear()
Definition: eval_points.hh:117
IntegralPtrMap< internal_integrals::Bulk > bulk_integrals_
Maps of all BulkIntegrals of dimensions 0,1,2,3.
Definition: eval_points.hh:204
uint get_max_bulk_quad_size(unsigned int dim) const
Return maximal size of quadrature of given dimension of bulk integral (Bulk, Coupling /lower-dim/)
Definition: eval_points.cc:89
int subset_size(unsigned int dim, unsigned int idx) const
Return number of local points corresponding to subset.
Definition: eval_points.hh:91
IntegralPtrMap< internal_integrals::Edge > edge_integrals_
Maps of all EdgeIntegrals of dimensions 1,2,3.
Definition: eval_points.hh:207
std::array< DimEvalPoints, 4 > dim_eval_points_
Sub objects of dimensions 0,1,2,3.
Definition: eval_points.hh:201
static constexpr unsigned int max_subsets
Maximal number of hold subsets.
Definition: eval_points.hh:60
arma::vec::fixed< dim > local_point(unsigned int local_point_idx) const
Return local coordinates of given local point and appropriate dim.
Definition: eval_points.hh:75
EvalPoints()
Constructor.
Definition: eval_points.cc:29
uint get_max_side_quad_size(unsigned int dim) const
Return maximal size of quadrature of given dimension of side integral (Edge, Coupling /higher-dim/,...
Definition: eval_points.cc:93
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
std::unordered_map< std::tuple< uint, uint >, std::shared_ptr< Integral >, IntegralTplHash > IntegralPtrMap
Alias for unordered_map of shared_ptr<Integral> with custom hash.
unsigned int uint
Implementation of range helper class.