Flow123d  jenkins-Flow123d-windows-release-multijob-285
mapping.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Class Mapping calculates data related to the mapping
27  * of the reference cell to the actual cell, such as Jacobian
28  * and normal vectors.
29  * @author Jan Stebel
30  */
31 
32 #ifndef MAPPING_HH_
33 #define MAPPING_HH_
34 
35 #include <armadillo>
36 #include <vector>
37 #include "fem/dofhandler.hh"
38 #include "mesh/side_impl.hh"
39 #include "fem/update_flags.hh"
40 #include "mesh/ref_element.hh"
41 
42 
43 
44 template<unsigned int dim> class Quadrature;
45 template<unsigned int dim, unsigned int spacedim> class FEValuesData;
46 
47 
48 
49 
50 /**
51  * @brief Calculates determinant of a rectangular matrix.
52  */
53 template<class T>
54 double determinant(const T &M);
55 
56 
57 template<> inline double determinant(const arma::mat::fixed<1,2> &M)
58 {
59  return sqrt(M(0,0)*M(0,0)+M(0,1)*M(0,1));
60 }
61 
62 template<> inline double determinant(const arma::mat::fixed<2,1> &M)
63 {
64  return sqrt(M(0,0)*M(0,0)+M(1,0)*M(1,0));
65 }
66 
67 template<> inline double determinant(const arma::mat::fixed<0,3> &M)
68 {
69  return 0;
70 }
71 
72 template<> inline double determinant(const arma::mat::fixed<3,0> &M)
73 {
74  return 0;
75 }
76 
77 template<> inline double determinant(const arma::mat::fixed<1,3> &M)
78 {
79  return sqrt(M(0,0)*M(0,0)+M(0,1)*M(0,1)+M(0,2)*M(0,2));
80 }
81 
82 template<> inline double determinant(const arma::mat::fixed<3,1> &M)
83 {
84  return sqrt(M(0,0)*M(0,0)+M(1,0)*M(1,0)+M(2,0)*M(2,0));
85 }
86 
87 template<> inline double determinant(const arma::mat::fixed<2,3> &M)
88 {
89  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))
90  -(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)));
91 }
92 
93 template<> inline double determinant(const arma::mat::fixed<3,2> &M)
94 {
95  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))
96  -(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)));
97 }
98 
99 template<unsigned int n> inline double determinant(const arma::mat::fixed<n,n> &M)
100 {
101  return det(M);
102 }
103 
104 
105 /**
106  * @brief Mapping data that can be precomputed on the actual cell.
107  *
108  * So far this involves only the (local) barycentric coordinates of quadrature points.
109  */
111 {
112 public:
113  /**
114  * @brief Auxiliary array of barycentric coordinates of quadrature points.
115  */
117 };
118 
119 
120 
121 /**
122  * @brief Abstract class for the mapping between reference and actual cell.
123  *
124  * Class Mapping calculates data related to the mapping of the
125  * reference cell to the actual cell, such as Jacobian and normal
126  * vectors.
127  */
128 template<unsigned int dim, unsigned int spacedim>
129 class Mapping
130 {
131 public:
132 
133  /**
134  * @brief Calculates the mapping data on the reference cell.
135  *
136  * @param q Quadrature rule.
137  * @param flags Update flags.
138  */
139  virtual MappingInternalData *initialize(const Quadrature<dim> &q, UpdateFlags flags) = 0;
140 
141  /**
142  * @brief Decides which additional quantities have to be computed
143  * for each cell.
144  *
145  * @param flags Flags of required quantities.
146  * @return Flags of all necessary quantities.
147  */
148  virtual UpdateFlags update_each(UpdateFlags flags) = 0;
149 
150  /**
151  * @brief Calculates the mapping data and stores them in the provided
152  * structures.
153  *
154  * @param cell The actual cell.
155  * @param q Quadrature rule.
156  * @param data Precomputed mapping data.
157  * @param fv_data Data to be computed.
158  */
159  virtual void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell,
160  const Quadrature<dim> &q,
161  MappingInternalData &data,
162  FEValuesData<dim,spacedim> &fv_data) = 0;
163 
164  /**
165  * @brief Calculates the mapping data related to a given side, namely the
166  * jacobian determinants and the normal vectors.
167  *
168  * @param cell The actual cell.
169  * @param side Number of the side.
170  * @param q Quadrature rule.
171  * @param data Precomputed mapping data.
172  * @param fv_data Data to be computed.
173  */
174  virtual void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &cell,
175  unsigned int sid,
176  const Quadrature<dim> &q,
177  MappingInternalData &data,
178  FEValuesData<dim,spacedim> &fv_data) = 0;
179 
180  /**
181  * @brief Creates a cell dim-dimensional quadrature from side (dim-1)-dimensional quadrature.
182  *
183  * @param cell The actual cell.
184  * @param sid The side id.
185  * @param pid Permutations index.
186  * @param subq @p dim-1 dimensional quadrature for integration on the side.
187  * @param q The computed @p dim dimensional quadrature.
188  */
189  void transform_subquadrature(unsigned int sid,
190  unsigned int pid,
191  const Quadrature<dim-1> &subq,
192  Quadrature<dim> &q);
193 
194  /// Destructor.
195  virtual ~Mapping() {};
196 };
197 
198 
199 
200 template<unsigned int dim, unsigned int spacedim> inline
202  unsigned int pid,
203  const Quadrature<dim - 1> & subq,
204  Quadrature<dim> &q)
205 {
206  q.resize(subq.size());
207 
208  double lambda;
209 
210  // vectors of barycentric coordinates of quadrature points
211  arma::vec::fixed<dim+1> el_bar_coords;
212  arma::vec::fixed<dim> side_bar_coords;
213 
214  for (unsigned int k=0; k<q.size(); k++)
215  {
216  // Calculate barycentric coordinates on the side of the k-th
217  // quadrature point.
218  el_bar_coords.zeros();
219  lambda = 0;
220  // Apply somewhere permutation of indices!
221  for (int j=0; j<dim-1; j++)
222  {
223  side_bar_coords(j) = (subq.point(k))[j];
224  lambda += (subq.point(k))[j];
225  }
226  side_bar_coords(dim-1) = 1.-lambda;
227 
228  // transform to element coordinates
229  for (unsigned int i=0; i<dim; i++)
230  el_bar_coords((RefElement<dim>::side_nodes[sid][RefElement<dim>::side_permutations[pid][i]]+dim)%(dim+1)) = side_bar_coords((i+dim-1)%dim);
231  q.set_point(k, el_bar_coords.subvec(0,dim-1));
232  q.set_weight(k, subq.weight(k));
233  }
234 }
235 
236 
237 
238 
239 
240 
241 
242 
243 
244 
245 
246 
247 #endif /* MAPPING_HH_ */
void set_weight(const unsigned int i, const double w)
Sets individual quadrature weight.
Definition: quadrature.hh:168
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:78
virtual void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &cell, unsigned int sid, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)=0
Calculates the mapping data related to a given side, namely the jacobian determinants and the normal ...
virtual ~Mapping()
Destructor.
Definition: mapping.hh:195
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
virtual MappingInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)=0
Calculates the mapping data on the reference cell.
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:42
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:158
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:44
void resize(const unsigned int n_q_points)
Modify the number of quadrature points.
Definition: quadrature.hh:127
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:201
std::vector< arma::vec > bar_coords
Auxiliary array of barycentric coordinates of quadrature points.
Definition: mapping.hh:116
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:57
void set_point(const unsigned int i, const arma::vec::fixed< dim > &p)
Sets individual quadrature point coordinates.
Definition: quadrature.hh:152
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:136
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:141
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
virtual void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)=0
Calculates the mapping data and stores them in the provided structures.
virtual UpdateFlags update_each(UpdateFlags flags)=0
Decides which additional quantities have to be computed for each cell.
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:110