Flow123d  release_3.0.0-885-g06276d1
quadrature.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 quadrature.hh
15  * @brief Basic definitions of numerical quadrature rules.
16  * @author Jan Stebel
17  */
18 
19 #ifndef QUADRATURE_HH_
20 #define QUADRATURE_HH_
21 
22 #include <armadillo>
23 #include <vector>
24 
25 #include "mesh/ref_element.hh"
26 
27 /**
28  * @brief Base class for quadrature rules on simplices in arbitrary dimensions.
29  *
30  * This class stores quadrature points and weights on the reference line,
31  * triangle or tetrahedron, respectively.
32  * Quadrature rules are used for evaluation of integrals over elements.
33  * In particular, for a reference element @f$E@f$ we have:
34  * @f[
35  * \int_E f(x)\,dx \approx \sum_{i=1}^{N} w_i f(p_i),
36  * @f]
37  * where @f$\{w_i\}@f$, @f$\{p_i\}@f$ are the quadrature weights and the quadrature points,
38  * respectively.
39  *
40  * TODO:
41  * - remove set_weight, set point; quadrature should be set by its descendants
42  * - introduce Quadrature point which should store point coords together with the weight in raw double[1+dim] array
43  * and just return arma object on the flys
44  */
45 template<unsigned int dim>
46 class Quadrature {
47 public:
48  /**
49  * @brief Constructor.
50  * @param n_quadrature_points Number of quadrature points to be allocated.
51  */
52  Quadrature(const unsigned int n_quadrature_points = 0);
53 
54  /// Copy constructor.
55  Quadrature(const Quadrature<dim> &q);
56 
57  /** @brief Constructor from quadrature of lower dimension (e.g. for side integration).
58  * @param sub_quadrature lower dimnesional (dim-1) quadrature
59  * @param sid local index of side
60  * @param pid index of permutation of nodes on given side
61  */
62  Quadrature(const Quadrature<dim-1> &sub_quadrature, unsigned int sid, unsigned int pid);
63 
64  /// Virtual destructor.
65  virtual ~Quadrature();
66 
67  /**
68  * @brief Modify the number of quadrature points.
69  * @param n_q_points New number of quadrature points.
70  */
71  void resize(const unsigned int n_q_points);
72 
73  /// Returns number of quadrature points.
74  const unsigned int size() const;
75 
76  /// Returns the <tt>i</tt>th quadrature point.
77  const arma::vec::fixed<dim> & point(const unsigned int i) const;
78 
79  /// Return a reference to the whole array of quadrature points.
81 
82  /**
83  * @brief Sets individual quadrature point coordinates.
84  * @param i Number of the quadrature point.
85  * @param p New coordinates.
86  */
87  void set_point(const unsigned int i, const arma::vec::fixed<dim> &p);
88 
89  /// Returns the <tt>i</tt>th weight.
90  double weight(const unsigned int i) const;
91 
92  /// Return a reference to the whole array of weights.
93  const std::vector<double> & get_weights() const;
94 
95  /// Sets individual quadrature weight.
96  void set_weight(const unsigned int i, const double w);
97 
98 protected:
99  /**
100  * @brief List of quadrature points.
101  *
102  * To be filled by the constructors of the derived classes.
103  */
105 
106  /**
107  * @brief List of weights to the quadrature points.
108  *
109  * To be filled by the constructors of the derived classes.
110  */
112 
113 };
114 
115 
116 
117 template<unsigned int dim>
118 Quadrature<dim>::Quadrature(const unsigned int n_q)
119 {
120  resize(n_q);
121 }
122 
123 template<unsigned int dim>
126  weights(q.weights)
127 {}
128 
129 template<unsigned int dim>
130 void Quadrature<dim>::resize(const unsigned int n_q)
131 {
132  arma::vec::fixed<dim> v;
133  v.fill(0);
134  quadrature_points.resize(n_q, v);
135  weights.resize(n_q, 0);
136 }
137 
138 template<unsigned int dim>
139 inline const unsigned int Quadrature<dim>::size() const {
140  return weights.size();
141 }
142 
143 template<unsigned int dim>
144 inline const arma::vec::fixed<dim> & Quadrature<dim>::point(
145  const unsigned int i) const {
146  return quadrature_points[i];
147 }
148 
149 template<unsigned int dim>
151  return quadrature_points;
152 }
153 
154 template<unsigned int dim>
155 inline void Quadrature<dim>::set_point(const unsigned int i, const arma::vec::fixed<dim> &p)
156 {
157  quadrature_points[i] = p;
158 }
159 
160 template<unsigned int dim>
161 inline double Quadrature<dim>::weight(const unsigned int i) const {
162  return weights[i];
163 }
164 
165 template<unsigned int dim>
167  return weights;
168 }
169 
170 template<unsigned int dim>
171 inline void Quadrature<dim>::set_weight(const unsigned int i, const double w)
172 {
173  weights[i] = w;
174 }
175 
176 template<unsigned int dim>
178 {}
179 
180 template<unsigned int dim> inline
181 Quadrature<dim>::Quadrature(const Quadrature<dim-1> &subq, unsigned int sid, unsigned int pid)
182 {
183  resize(subq.size());
184 
185 // double lambda;
186 //
187 // // vectors of barycentric coordinates of quadrature points
188 // arma::vec::fixed<dim+1> el_bar_coords;
189 // arma::vec::fixed<dim> side_bar_coords;
190 //
191 // for (unsigned int k=0; k<subq.size(); k++)
192 // {
193 // const arma::vec::fixed<dim-1> &sub_point = subq.point(k);
194 // // Calculate barycentric coordinates on the side of the k-th
195 // // quadrature point.
196 // el_bar_coords.zeros();
197 // lambda = 0;
198 // // Apply somewhere permutation of indices!
199 // for (unsigned int j=0; j<dim-1; j++)
200 // {
201 // side_bar_coords(j) = sub_point(j);
202 // lambda += sub_point(j);
203 // }
204 // side_bar_coords(dim-1) = 1.0 - lambda;
205 //
206 // // transform to element coordinates
207 // auto side_nodes = RefElement<dim>::interact(Interaction<0, (dim - 1)>(sid));
208 // for (unsigned int i=0; i<dim; i++) {
209 // // TODO: use RefElement<>::interpolate to map coordinates from the subelement
210 // unsigned int i_node = (side_nodes[RefElement<dim>::side_permutations[pid][i]]+dim)%(dim+1);
211 // el_bar_coords(i_node) = side_bar_coords((i+dim-1)%dim);
212 // }
213 // quadrature_points[k] = el_bar_coords.subvec(0,dim-1);
214 // weights[k] = subq.weight(k);
215 // }
216 
217  arma::vec::fixed<dim+1> el_bar_coords, final_bar;
218 
219  for (unsigned int k=0; k<subq.size(); k++)
220  {
221  //compute barycentric coordinates on element
222  arma::vec::fixed<dim> p = RefElement<dim-1>::local_to_bary(subq.point(k));
223  arma::vec::fixed<dim> pp;
224 
225  //permute
226  for (unsigned int i=0; i<RefElement<dim>::n_nodes_per_side; i++) {
227  pp(i) = p(RefElement<dim>::side_permutations[pid][i]);
228  }
229 
230  el_bar_coords = RefElement<dim>::template interpolate<dim-1>(pp,sid);
231 
232  //get local coordinates and set
234  weights[k] = subq.weight(k);
235  }
236 }
237 
238 // Specialized subquadrature consructor for dim=1.
239 template<> inline
240 Quadrature<1>::Quadrature(const Quadrature<0> &subq, unsigned int sid, unsigned int pid)
241 {
242  arma::vec::fixed<1> p({(double)sid});
243  quadrature_points.push_back(p);
244  weights.push_back(1);
245 }
246 
247 #endif /* QUADRATURE_HH_ */
void set_weight(const unsigned int i, const double w)
Sets individual quadrature weight.
Definition: quadrature.hh:171
std::vector< double > weights
List of weights to the quadrature points.
Definition: quadrature.hh:111
Quadrature(const unsigned int n_quadrature_points=0)
Constructor.
Definition: quadrature.hh:118
std::vector< arma::vec::fixed< dim > > quadrature_points
List of quadrature points.
Definition: quadrature.hh:104
const std::vector< arma::vec::fixed< dim > > & get_points() const
Return a reference to the whole array of quadrature points.
Definition: quadrature.hh:150
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:35
static LocalPoint bary_to_local(const BaryPoint &bp)
Converts from barycentric to local coordinates.
Definition: ref_element.cc:300
virtual ~Quadrature()
Virtual destructor.
Definition: quadrature.hh:177
double weight(const unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:161
void resize(const unsigned int n_q_points)
Modify the number of quadrature points.
Definition: quadrature.hh:130
void set_point(const unsigned int i, const arma::vec::fixed< dim > &p)
Sets individual quadrature point coordinates.
Definition: quadrature.hh:155
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:139
const std::vector< double > & get_weights() const
Return a reference to the whole array of weights.
Definition: quadrature.hh:166
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:144
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
Definition: ref_element.cc:286
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...