Flow123d  release_3.0.0-1263-g7cf53c1
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/accessors.hh"
28 #include "fem/update_flags.hh"
29 #include "mesh/ref_element.hh"
30 #include "system/fmt/posix.h" // for FMT_UNUSED
31 
32 
33 
34 class Quadrature;
35 template<unsigned int dim, unsigned int spacedim> class FEValuesData;
36 
37 
38 
39 
40 /**
41  * @brief Calculates determinant of a rectangular matrix.
42  */
43 template<class T>
44 double determinant(const T &M);
45 
46 
47 template<> inline double determinant(const arma::mat::fixed<1,2> &M)
48 {
49  return sqrt(M(0,0)*M(0,0)+M(0,1)*M(0,1));
50 }
51 
52 template<> inline double determinant(const arma::mat::fixed<2,1> &M)
53 {
54  return sqrt(M(0,0)*M(0,0)+M(1,0)*M(1,0));
55 }
56 
57 template<> inline double determinant(FMT_UNUSED const arma::mat::fixed<0,3> &M)
58 {
59  return 0;
60 }
61 
62 template<> inline double determinant(FMT_UNUSED const arma::mat::fixed<3,0> &M)
63 {
64  return 0;
65 }
66 
67 template<> inline double determinant(const arma::mat::fixed<1,3> &M)
68 {
69  return sqrt(M(0,0)*M(0,0)+M(0,1)*M(0,1)+M(0,2)*M(0,2));
70 }
71 
72 template<> inline double determinant(const arma::mat::fixed<3,1> &M)
73 {
74  return sqrt(M(0,0)*M(0,0)+M(1,0)*M(1,0)+M(2,0)*M(2,0));
75 }
76 
77 template<> inline double determinant(const arma::mat::fixed<2,3> &M)
78 {
79  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))
80  -(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)));
81 }
82 
83 template<> inline double determinant(const arma::mat::fixed<3,2> &M)
84 {
85  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))
86  -(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)));
87 }
88 
89 template<arma::uword n> inline double determinant(const arma::mat::fixed<n,n> &M)
90 {
91  return det(M);
92 }
93 
94 
95 /**
96  * @brief Mapping data that can be precomputed on the actual cell.
97  *
98  * So far this involves only the (local) barycentric coordinates of quadrature points.
99  */
101 {
102 public:
103  /**
104  * @brief Auxiliary array of barycentric coordinates of quadrature points.
105  */
107 };
108 
109 
110 
111 /**
112  * @brief Abstract class for the mapping between reference and actual cell.
113  *
114  * Class Mapping calculates data related to the mapping of the
115  * reference cell to the actual cell, such as Jacobian and normal
116  * vectors.
117  */
118 template<unsigned int dim, unsigned int spacedim>
119 class Mapping
120 {
121 public:
122 
123  /**
124  * @brief Calculates the mapping data on the reference cell.
125  *
126  * @param q Quadrature rule.
127  * @param flags Update flags.
128  */
129  virtual MappingInternalData *initialize(const Quadrature &q, UpdateFlags flags) = 0;
130 
131  /**
132  * @brief Decides which additional quantities have to be computed
133  * for each cell.
134  *
135  * @param flags Flags of required quantities.
136  * @return Flags of all necessary quantities.
137  */
138  virtual UpdateFlags update_each(UpdateFlags flags) = 0;
139 
140  /**
141  * @brief Calculates the mapping data and stores them in the provided
142  * structures.
143  *
144  * @param cell The actual cell.
145  * @param q Quadrature rule.
146  * @param data Precomputed mapping data.
147  * @param fv_data Data to be computed.
148  */
149  virtual void fill_fe_values(const ElementAccessor<3> &cell,
150  const Quadrature &q,
151  MappingInternalData &data,
152  FEValuesData<dim,spacedim> &fv_data) = 0;
153 
154  /**
155  * @brief Calculates the mapping data related to a given side, namely the
156  * jacobian determinants and the normal vectors.
157  *
158  * @param cell The actual cell.
159  * @param side Number of the side.
160  * @param q Quadrature rule.
161  * @param data Precomputed mapping data.
162  * @param fv_data Data to be computed.
163  */
164  virtual void fill_fe_side_values(const ElementAccessor<3> &cell,
165  unsigned int sid,
166  const Quadrature &q,
167  MappingInternalData &data,
168  FEValuesData<dim,spacedim> &fv_data) = 0;
169 
170  /// Destructor.
171  virtual ~Mapping() {};
172 };
173 
174 #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:171
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: quadrature.hh:48
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
#define FMT_UNUSED
Definition: posix.h:75
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:106
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:100