Flow123d  release_3.0.0-684-g928e266
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<arma::uword 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 ElementAccessor<3> &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 ElementAccessor<3> &cell,
164  unsigned int sid,
165  const Quadrature<dim> &q,
166  MappingInternalData &data,
167  FEValuesData<dim,spacedim> &fv_data) = 0;
168 
169  /// Destructor.
170  virtual ~Mapping() {};
171 };
172 
173 #endif /* MAPPING_HH_ */
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:170
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:35
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:38
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:86
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