Flow123d  release_2.2.1-10-gb9fad4d
mapping.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 mapping.hh
15  * @brief Class Mapping calculates data related to the mapping
16  * of the reference cell to the actual cell, such as Jacobian
17  * and normal vectors.
18  * @author Jan Stebel
19  */
20 
21 #ifndef MAPPING_HH_
22 #define MAPPING_HH_
23 
24 #include <armadillo>
25 #include <vector>
26 #include "fem/dofhandler.hh"
27 #include "mesh/side_impl.hh"
28 #include "fem/update_flags.hh"
29 #include "mesh/ref_element.hh"
30 
31 
32 
33 template<unsigned int dim> class Quadrature;
34 template<unsigned int dim, unsigned int spacedim> class FEValuesData;
35 
36 
37 
38 
39 /**
40  * @brief Calculates determinant of a rectangular matrix.
41  */
42 template<class T>
43 double determinant(const T &M);
44 
45 
46 template<> inline double determinant(const arma::mat::fixed<1,2> &M)
47 {
48  return sqrt(M(0,0)*M(0,0)+M(0,1)*M(0,1));
49 }
50 
51 template<> inline double determinant(const arma::mat::fixed<2,1> &M)
52 {
53  return sqrt(M(0,0)*M(0,0)+M(1,0)*M(1,0));
54 }
55 
56 template<> inline double determinant(const arma::mat::fixed<0,3> &M)
57 {
58  return 0;
59 }
60 
61 template<> inline double determinant(const arma::mat::fixed<3,0> &M)
62 {
63  return 0;
64 }
65 
66 template<> inline double determinant(const arma::mat::fixed<1,3> &M)
67 {
68  return sqrt(M(0,0)*M(0,0)+M(0,1)*M(0,1)+M(0,2)*M(0,2));
69 }
70 
71 template<> inline double determinant(const arma::mat::fixed<3,1> &M)
72 {
73  return sqrt(M(0,0)*M(0,0)+M(1,0)*M(1,0)+M(2,0)*M(2,0));
74 }
75 
76 template<> inline double determinant(const arma::mat::fixed<2,3> &M)
77 {
78  return sqrt((M(0,0)*M(0,0)+M(0,1)*M(0,1)+M(0,2)*M(0,2))*(M(1,0)*M(1,0)+M(1,1)*M(1,1)+M(1,2)*M(1,2))
79  -(M(0,0)*M(1,0)+M(0,1)*M(1,1)+M(0,2)*M(1,2))*(M(0,0)*M(1,0)+M(0,1)*M(1,1)+M(0,2)*M(1,2)));
80 }
81 
82 template<> inline double determinant(const arma::mat::fixed<3,2> &M)
83 {
84  return sqrt((M(0,0)*M(0,0)+M(1,0)*M(1,0)+M(2,0)*M(2,0))*(M(0,1)*M(0,1)+M(1,1)*M(1,1)+M(2,1)*M(2,1))
85  -(M(0,0)*M(0,1)+M(1,0)*M(1,1)+M(2,0)*M(2,1))*(M(0,0)*M(0,1)+M(1,0)*M(1,1)+M(2,0)*M(2,1)));
86 }
87 
88 template<unsigned int n> inline double determinant(const arma::mat::fixed<n,n> &M)
89 {
90  return det(M);
91 }
92 
93 
94 /**
95  * @brief Mapping data that can be precomputed on the actual cell.
96  *
97  * So far this involves only the (local) barycentric coordinates of quadrature points.
98  */
100 {
101 public:
102  /**
103  * @brief Auxiliary array of barycentric coordinates of quadrature points.
104  */
106 };
107 
108 
109 
110 /**
111  * @brief Abstract class for the mapping between reference and actual cell.
112  *
113  * Class Mapping calculates data related to the mapping of the
114  * reference cell to the actual cell, such as Jacobian and normal
115  * vectors.
116  */
117 template<unsigned int dim, unsigned int spacedim>
118 class Mapping
119 {
120 public:
121 
122  /**
123  * @brief Calculates the mapping data on the reference cell.
124  *
125  * @param q Quadrature rule.
126  * @param flags Update flags.
127  */
128  virtual MappingInternalData *initialize(const Quadrature<dim> &q, UpdateFlags flags) = 0;
129 
130  /**
131  * @brief Decides which additional quantities have to be computed
132  * for each cell.
133  *
134  * @param flags Flags of required quantities.
135  * @return Flags of all necessary quantities.
136  */
137  virtual UpdateFlags update_each(UpdateFlags flags) = 0;
138 
139  /**
140  * @brief Calculates the mapping data and stores them in the provided
141  * structures.
142  *
143  * @param cell The actual cell.
144  * @param q Quadrature rule.
145  * @param data Precomputed mapping data.
146  * @param fv_data Data to be computed.
147  */
148  virtual void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell,
149  const Quadrature<dim> &q,
150  MappingInternalData &data,
151  FEValuesData<dim,spacedim> &fv_data) = 0;
152 
153  /**
154  * @brief Calculates the mapping data related to a given side, namely the
155  * jacobian determinants and the normal vectors.
156  *
157  * @param cell The actual cell.
158  * @param side Number of the side.
159  * @param q Quadrature rule.
160  * @param data Precomputed mapping data.
161  * @param fv_data Data to be computed.
162  */
163  virtual void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &cell,
164  unsigned int sid,
165  const Quadrature<dim> &q,
166  MappingInternalData &data,
167  FEValuesData<dim,spacedim> &fv_data) = 0;
168 
169  /**
170  * @brief Creates a cell dim-dimensional quadrature from side (dim-1)-dimensional quadrature.
171  *
172  * @param cell The actual cell.
173  * @param sid The side id.
174  * @param pid Permutations index.
175  * @param subq @p dim-1 dimensional quadrature for integration on the side.
176  * @param q The computed @p dim dimensional quadrature.
177  */
178  void transform_subquadrature(unsigned int sid,
179  unsigned int pid,
180  const Quadrature<dim-1> &subq,
181  Quadrature<dim> &q);
182 
183  /// Destructor.
184  virtual ~Mapping() {};
185 };
186 
187 
188 
189 template<unsigned int dim, unsigned int spacedim> inline
191  unsigned int pid,
192  const Quadrature<dim - 1> & subq,
193  Quadrature<dim> &q)
194 {
195  q.resize(subq.size());
196 
197  double lambda;
198 
199  // vectors of barycentric coordinates of quadrature points
200  arma::vec::fixed<dim+1> el_bar_coords;
201  arma::vec::fixed<dim> side_bar_coords;
202 
203  for (unsigned int k=0; k<q.size(); k++)
204  {
205  // Calculate barycentric coordinates on the side of the k-th
206  // quadrature point.
207  el_bar_coords.zeros();
208  lambda = 0;
209  // Apply somewhere permutation of indices!
210  for (int j=0; j<dim-1; j++)
211  {
212  side_bar_coords(j) = (subq.point(k))[j];
213  lambda += (subq.point(k))[j];
214  }
215  side_bar_coords(dim-1) = 1.-lambda;
216 
217  // transform to element coordinates
218  for (unsigned int i=0; i<dim; i++)
219  el_bar_coords((RefElement<dim>::side_nodes[sid][RefElement<dim>::side_permutations[pid][i]]+dim)%(dim+1)) = side_bar_coords((i+dim-1)%dim);
220  q.set_point(k, el_bar_coords.subvec(0,dim-1));
221  q.set_weight(k, subq.weight(k));
222  }
223 }
224 
225 
226 
227 
228 
229 
230 
231 
232 
233 
234 
235 
236 #endif /* MAPPING_HH_ */
void set_weight(const unsigned int i, const double w)
Sets individual quadrature weight.
Definition: quadrature.hh:157
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
virtual ~Mapping()
Destructor.
Definition: mapping.hh:184
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:31
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
double weight(const unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:147
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:33
void resize(const unsigned int n_q_points)
Modify the number of quadrature points.
Definition: quadrature.hh:116
void transform_subquadrature(unsigned int sid, unsigned int pid, const Quadrature< dim-1 > &subq, Quadrature< dim > &q)
Creates a cell dim-dimensional quadrature from side (dim-1)-dimensional quadrature.
Definition: mapping.hh:190
std::vector< arma::vec > bar_coords
Auxiliary array of barycentric coordinates of quadrature points.
Definition: mapping.hh:105
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:46
void set_point(const unsigned int i, const arma::vec::fixed< dim > &p)
Sets individual quadrature point coordinates.
Definition: quadrature.hh:141
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:125
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:130
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:99