Flow123d  release_3.0.0-1193-g9220a69
mapping_p1.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_p1.hh
15  * @brief Class MappingP1 implements the affine transformation of
16  * the unit cell onto the actual cell.
17  * @author Jan Stebel
18  */
19 
20 #ifndef MAPPING_P1_HH_
21 #define MAPPING_P1_HH_
22 
23 #include <string.h> // for memcpy
24 #include <algorithm> // for min, swap
25 #include <cmath> // for abs, pow, fabs, log10
26 #include <ostream> // for operator<<
27 #include <string> // for operator<<
28 #include <armadillo>
29 #include "fem/dofhandler.hh" // for DOFHandlerBase, DOFHan...
30 #include "fem/mapping.hh" // for MappingInternalData (p...
31 #include "fem/update_flags.hh" // for operator&, operator|
32 #include "mesh/elements.h" // for Element::side
33 #include "mesh/side_impl.hh" // for Side::node
34 template <unsigned int dim, unsigned int spacedim> class FEValuesData;
35 class Quadrature;
36 
37 
38 /**
39  * Auxiliary templates to resolve nonzero matrix sizes for small dimensions.
40  * ( some compilers do not accept ternary operator in template parameters.)
41  */
42 template<unsigned int dim>
43 class MatrixSizes {
44 public:
45  static const unsigned int dim_minus_one = dim-1;
46 };
47 
48 template<>
49 class MatrixSizes<0> {
50 public:
51  static const unsigned int dim_minus_one = 0;
52 };
53 
54 
55 /**
56  * @brief Affine mapping between reference and actual cell.
57  *
58  * Class MappingP1 implements the affine transformation of
59  * the reference cell onto the actual cell.
60  *
61  * @param dim Dimension of the cells.
62  * @param spacedim Dimension of the Euclidean space.
63  */
64 template<unsigned int dim, unsigned int spacedim>
65 class MappingP1 : public Mapping<dim,spacedim>
66 {
67 public:
68 
69  typedef arma::vec::fixed<dim+1> BaryPoint;
70  typedef arma::vec::fixed<spacedim> RealPoint;
71  typedef arma::mat::fixed<spacedim, dim+1> ElementMap;
72 
73  /**
74  * @brief Constructor.
75  */
76  MappingP1();
77 
78  /**
79  * @brief Initializes the structures and computes static data.
80  *
81  * @param q Quadrature rule.
82  * @param flags Update flags.
83  * @return The computed mapping data.
84  */
85  MappingInternalData *initialize(const Quadrature &q, UpdateFlags flags);
86 
87  /**
88  * @brief Determines which additional quantities have to be computed.
89  *
90  * @param flags Update flags for required quantities.
91  * @return All necessary flags.
92  */
93  UpdateFlags update_each(UpdateFlags flags);
94 
95  /**
96  * @brief Calculates the mapping data on the actual cell.
97  *
98  * @param cell The actual cell.
99  * @param q Quadrature rule.
100  * @param data Precomputed mapping data.
101  * @param fv_data Data to be computed.
102  */
103  void fill_fe_values(const ElementAccessor<3> &cell,
104  const Quadrature &q,
105  MappingInternalData &data,
106  FEValuesData<dim,spacedim> &fv_data);
107 
108  /**
109  * @brief Calculates the mapping data on a side of a cell.
110  *
111  * @param cell The actual cell.
112  * @param sid Number of the side.
113  * @param q The quadrature rule with points on the side.
114  * @param data Precomputed mapping data.
115  * @param fv_data Data to be computed.
116  */
117  void fill_fe_side_values(const ElementAccessor<3> &cell,
118  unsigned int sid,
119  const Quadrature &q,
120  MappingInternalData &data,
121  FEValuesData<dim,spacedim> &fv_data);
122 
123 
124  /**
125  * Map from reference element (barycentric coords) to global coord system.
126  * Matrix(3, dim+1) M: x_real = M * x_bary;
127  * M columns are real coordinates of nodes.
128  */
129  ElementMap element_map(ElementAccessor<3> elm) const;
130 
131  /**
132  * Project given point in real coordinates to reference element (barycentic coordinates).
133  * Result vector have dimension dim()+1.
134  * Use RefElement<dim>::bary_to_local() to get local coordinates.
135  */
136  BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map) const;
137 
138  /**
139  * Project given point from reference element (barycentic coordinates) to real coordinates.
140  * Use RefElement<dim>::local_to_bary() to get barycentric coordinates in input.
141  */
142  RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map) const;
143 
144  /**
145  * Clip a point given by barycentric cocordinates to the element.
146  * If the point is out of the element the closest point
147  * projection to the element surface is used.
148  */
149  BaryPoint clip_to_element(BaryPoint &barycentric);
150 
151  /// Test if element contains given point.
152  bool contains_point(arma::vec point, ElementAccessor<3> elm);
153 
154 private:
155 
156  /**
157  * @brief Auxiliary matrix of gradients of shape functions (used for
158  * computation of the Jacobian).
159  */
160  arma::mat::fixed<dim+1,dim> grad;
161 
162 };
163 
164 
165 
166 
167 
168 
169 
170 
171 #endif /* MAPPING_P1_HH_ */
arma::mat::fixed< spacedim, dim+1 > ElementMap
Definition: mapping_p1.hh:71
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
static const unsigned int dim_minus_one
Definition: mapping_p1.hh:45
Mat< double, N, 1 > vec
Definition: armor.hh:211
arma::vec::fixed< spacedim > RealPoint
Definition: mapping_p1.hh:70
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
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:38
arma::vec::fixed< dim+1 > BaryPoint
Definition: mapping_p1.hh:69
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:65
Class Mapping calculates data related to the mapping of the reference cell to the actual cell...
arma::mat::fixed< dim+1, dim > grad
Auxiliary matrix of gradients of shape functions (used for computation of the Jacobian).
Definition: mapping_p1.hh:160
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:86
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:99