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