Flow123d  JS_before_hm-1730-gbf8dc60dc
eval_points.cc
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.cc
15  * @brief
16  * @author David Flanderka
17  */
18 
19 #include "fields/eval_points.hh"
20 #include "fields/eval_subset.hh"
21 #include "quadrature/quadrature.hh"
22 #include "mesh/ref_element.hh"
23 #include <memory>
24 
25 
26 const unsigned int EvalPoints::undefined_dim = 10;
27 
29 : dim_eval_points_({DimEvalPoints(0), DimEvalPoints(1), DimEvalPoints(2), DimEvalPoints(3)}), max_size_(0)
30 {}
31 
32 template <unsigned int dim>
33 std::shared_ptr<BulkIntegral> EvalPoints::add_bulk(const Quadrature &quad)
34 {
35  ASSERT_EQ(dim, quad.dim());
36  std::shared_ptr<BulkIntegral> bulk_integral = std::make_shared<BulkIntegral>(shared_from_this(), dim);
37  dim_eval_points_[dim].add_local_points<dim>( quad.get_points() );
38  dim_eval_points_[dim].add_subset();
39  this->set_max_size();
40  return bulk_integral;
41 }
42 
43 template <>
44 std::shared_ptr<BulkIntegral> EvalPoints::add_bulk<0>(const Quadrature &quad)
45 {
46  ASSERT_EQ(0, quad.dim());
47  std::shared_ptr<BulkIntegral> bulk_integral = std::make_shared<BulkIntegral>(shared_from_this(), 0);
48  dim_eval_points_[0].add_subset();
49  this->set_max_size();
50  return bulk_integral;
51 }
52 
53 template <unsigned int dim>
54 std::shared_ptr<EdgeIntegral> EvalPoints::add_edge(const Quadrature &quad)
55 {
56  ASSERT_EQ(dim, quad.dim()+1);
57  unsigned int old_data_size=this->size(dim), new_data_size; // interval of side subset data
58  unsigned int points_per_side = quad.make_from_side<dim>(0, 0).get_points().size();
59  unsigned int n_side_permutations = RefElement<dim>::n_side_permutations;
60 
61  std::shared_ptr<EdgeIntegral> edge_integral = std::make_shared<EdgeIntegral>(shared_from_this(), dim, n_side_permutations, points_per_side);
62  unsigned int*** perm_indices = edge_integral->perm_indices_;
63 
64  // permutation 0
65  for (unsigned int i=0; i<dim+1; ++i) { // sides
66  Quadrature high_dim_q = quad.make_from_side<dim>(i, 0);
67  dim_eval_points_[dim].add_local_points<dim>( high_dim_q.get_points() );
68  }
69  dim_eval_points_[dim].add_subset();
70  new_data_size = this->size(dim);
71  unsigned int i_data=old_data_size;
72  for (unsigned int i_side=0; i_side<dim+1; ++i_side) {
73  for (unsigned int i_point=0; i_point<points_per_side; ++i_point) {
74  perm_indices[i_side][0][i_point] = i_data;
75  ++i_data;
76  }
77  }
78 
79  // permutation 1...N
80  for (unsigned int i_perm=1; i_perm<n_side_permutations; ++i_perm) {
81  for (unsigned int i_side=0; i_side<dim+1; ++i_side) {
82  Quadrature high_dim_q = quad.make_from_side<dim>(i_side, i_perm);
83  const Armor::Array<double> & quad_points = high_dim_q.get_points();
84  for (unsigned int i_point=0; i_point<quad_points.size(); ++i_point) {
85  perm_indices[i_side][i_perm][i_point] = dim_eval_points_[dim].find_permute_point<dim>( quad_points.vec<dim>(i_point), old_data_size, new_data_size );
86  }
87  }
88  }
89 
90  this->set_max_size();
91  return edge_integral;
92 }
93 
94 template <unsigned int dim>
95 std::shared_ptr<CouplingIntegral> EvalPoints::add_coupling(const Quadrature &quad) {
96  ASSERT_EQ(dim, quad.dim()+1);
97 
98  std::shared_ptr<BulkIntegral> bulk_integral = this->add_bulk<dim-1>(quad);
99  std::shared_ptr<EdgeIntegral> edge_integral = this->add_edge<dim>(quad);
100  return std::make_shared<CouplingIntegral>(edge_integral, bulk_integral);
101 }
102 
103 template <unsigned int dim>
104 std::shared_ptr<BoundaryIntegral> EvalPoints::add_boundary(const Quadrature &quad) {
105  ASSERT_EQ(dim, quad.dim()+1);
106  std::shared_ptr<BulkIntegral> bulk_integral = this->add_bulk<dim-1>(quad);
107  std::shared_ptr<EdgeIntegral> edge_integral = this->add_edge<dim>(quad);
108  return std::make_shared<BoundaryIntegral>(edge_integral, bulk_integral);
109 }
110 
112 : local_points_(dim), n_subsets_(0), dim_(dim)
113 {
114  subset_starts_[0] = 0;
116 }
117 
118 template <unsigned int dim>
120  ASSERT_GT(dim, 0).error("Dimension 0 not supported!\n");
121  unsigned int local_points_old_size = local_points_.size();
122  local_points_.resize(quad_points.size() + local_points_old_size);
123  for (unsigned int i=0; i<quad_points.size(); ++i) {
124  local_points_.set(i+local_points_old_size) = quad_points.vec<dim>(i);
125  }
126 }
127 
128 template <unsigned int dim>
129 unsigned int EvalPoints::DimEvalPoints::find_permute_point(arma::vec coords, unsigned int data_begin, unsigned int data_end) {
130  ASSERT_GT(dim, 0).error("Dimension 0 not supported!\n");
131  for (unsigned int loc_idx=data_begin; loc_idx<data_end; ++loc_idx) {
132  // Check if point exists in local points vector.
133  if ( arma::norm(coords-local_points_.vec<dim>(loc_idx), 2) < 4*std::numeric_limits<double>::epsilon() ) return loc_idx;
134  }
135 
136  ASSERT(false);
137  return 0;
138 }
139 
141  ASSERT_LT_DBG(n_subsets_, EvalPoints::max_subsets).error("Maximal number of subsets exceeded!\n");
142 
143  n_subsets_++;
144  subset_starts_[n_subsets_] = this->size();
145 }
146 
147 
148 template std::shared_ptr<BulkIntegral> EvalPoints::add_bulk<0>(const Quadrature &);
149 template std::shared_ptr<BulkIntegral> EvalPoints::add_bulk<1>(const Quadrature &);
150 template std::shared_ptr<BulkIntegral> EvalPoints::add_bulk<2>(const Quadrature &);
151 template std::shared_ptr<BulkIntegral> EvalPoints::add_bulk<3>(const Quadrature &);
152 template std::shared_ptr<EdgeIntegral> EvalPoints::add_edge<1>(const Quadrature &);
153 template std::shared_ptr<EdgeIntegral> EvalPoints::add_edge<2>(const Quadrature &);
154 template std::shared_ptr<EdgeIntegral> EvalPoints::add_edge<3>(const Quadrature &);
155 template std::shared_ptr<CouplingIntegral> EvalPoints::add_coupling<2>(const Quadrature &);
156 template std::shared_ptr<CouplingIntegral> EvalPoints::add_coupling<3>(const Quadrature &);
157 template std::shared_ptr<BoundaryIntegral> EvalPoints::add_boundary<1>(const Quadrature &);
158 template std::shared_ptr<BoundaryIntegral> EvalPoints::add_boundary<2>(const Quadrature &);
159 template std::shared_ptr<BoundaryIntegral> EvalPoints::add_boundary<3>(const Quadrature &);
160 template void EvalPoints::DimEvalPoints::add_local_points<1>(const Armor::Array<double> &);
161 template void EvalPoints::DimEvalPoints::add_local_points<2>(const Armor::Array<double> &);
162 template void EvalPoints::DimEvalPoints::add_local_points<3>(const Armor::Array<double> &);
163 template unsigned int EvalPoints::DimEvalPoints::find_permute_point<1>(arma::vec, unsigned int, unsigned int);
164 template unsigned int EvalPoints::DimEvalPoints::find_permute_point<2>(arma::vec, unsigned int, unsigned int);
165 template unsigned int EvalPoints::DimEvalPoints::find_permute_point<3>(arma::vec, unsigned int, unsigned int);
ref_element.hh
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
eval_subset.hh
ASSERT_GT
#define ASSERT_GT(a, b)
Definition of comparative assert macro (Greater Than)
Definition: asserts.hh:312
Armor::vec
ArmaVec< double, N > vec
Definition: armor.hh:885
EvalPoints::DimEvalPoints::add_local_points
void add_local_points(const Armor::Array< double > &quad_points)
Adds set of local point to local_points_ (bulk or side of given permutation).
Definition: eval_points.cc:119
RefElement
Definition: ref_element.hh:163
EvalPoints::set_max_size
void set_max_size()
Definition: eval_points.hh:173
EvalPoints::DimEvalPoints::subset_starts_
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:168
ASSERT
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:347
Armor::Array::vec
ArmaVec< Type, nr > vec(uint mat_index) const
Definition: armor.hh:821
eval_points.hh
EvalPoints::undefined_dim
static const unsigned int undefined_dim
Undefined dimension of new (empty) object.
Definition: eval_points.hh:46
Quadrature::size
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:87
EvalPoints::DimEvalPoints
Subobject holds evaluation points data of one dimension (0,1,2,3)
Definition: eval_points.hh:115
EvalPoints::EvalPoints
EvalPoints()
Constructor.
Definition: eval_points.cc:28
EvalPoints::max_subset_points
static const unsigned int max_subset_points
Maximal average number of points hold in subset.
Definition: eval_points.hh:52
EvalPoints::add_coupling
std::shared_ptr< CouplingIntegral > add_coupling(const Quadrature &)
The same as add_bulk but for points between side points of element of dim and bulk points of element ...
Definition: eval_points.cc:95
EvalPoints::add_boundary
std::shared_ptr< BoundaryIntegral > add_boundary(const Quadrature &)
The same as add_bulk but for edge points on boundary sides.
Definition: eval_points.cc:104
EvalPoints::add_edge
std::shared_ptr< EdgeIntegral > add_edge(const Quadrature &)
The same as add_bulk but for edge points on sides.
Definition: eval_points.cc:54
ASSERT_EQ
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
Quadrature::dim
unsigned int dim() const
Definition: quadrature.hh:73
Quadrature::get_points
const Armor::Array< double > & get_points() const
Return a reference to the whole array of quadrature points.
Definition: quadrature.hh:104
EvalPoints::dim_eval_points_
std::array< DimEvalPoints, 4 > dim_eval_points_
Sub objects of dimensions 0,1,2,3.
Definition: eval_points.hh:178
Quadrature::make_from_side
Quadrature make_from_side(unsigned int sid, unsigned int pid) const
Definition: quadrature.cc:47
Armor::Array::size
unsigned int size() const
Definition: armor.hh:728
EvalPoints::DimEvalPoints::add_subset
void add_subset()
Adds new subset and its end size to subset_starts_ array.
Definition: eval_points.cc:140
Armor::Array::reinit
void reinit(uint size)
Definition: armor.hh:698
EvalPoints::max_subsets
static constexpr unsigned int max_subsets
Maximal number of hold subsets.
Definition: eval_points.hh:49
quadrature.hh
Basic definitions of numerical quadrature rules.
ASSERT_LT_DBG
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
Armor::Array< double >
EvalPoints::DimEvalPoints::find_permute_point
unsigned int find_permute_point(arma::vec coords, unsigned int data_begin, unsigned int data_end)
Find position of local point (coords) in subvector of local points given by limits <data_begin,...
Definition: eval_points.cc:129
Quadrature
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
EvalPoints::DimEvalPoints::local_points_
Armor::Array< double > local_points_
Local coords of points vector.
Definition: eval_points.hh:167
EvalPoints::size
unsigned int size(unsigned int dim) const
Return size of evaluation points object (number of points).
Definition: eval_points.hh:58
EvalPoints::DimEvalPoints::DimEvalPoints
DimEvalPoints(unsigned int dim)
Constructor.
Definition: eval_points.cc:111
EvalPoints::add_bulk
std::shared_ptr< BulkIntegral > add_bulk(const Quadrature &)
Definition: eval_points.cc:33