Flow123d  DF_mechanic_bench-4968b1b
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 <armadillo>
24 #include "fem/fem_tools.hh" // for MappingInternalData (p...
25 #include "fem/update_flags.hh" // for operator&, operator|
26 #include "mesh/accessors.hh" // for ElementAccessor
27 
28 class Quadrature;
29 
30 
31 /**
32  * Auxiliary templates to resolve nonzero matrix sizes for small dimensions.
33  * ( some compilers do not accept ternary operator in template parameters.)
34  */
35 template<unsigned int dim>
36 class MatrixSizes {
37 public:
38  static const unsigned int dim_minus_one = dim-1;
39 };
40 
41 template<>
42 class MatrixSizes<0> {
43 public:
44  static const unsigned int dim_minus_one = 0;
45 };
46 
47 
48 /**
49  * @brief Affine mapping between reference and actual cell.
50  *
51  * Class MappingP1 implements the affine transformation of
52  * the reference cell onto the actual cell.
53  *
54  * @param dim Dimension of the cells.
55  * @param spacedim Dimension of the Euclidean space.
56  */
57 template<unsigned int dim, unsigned int spacedim = 3>
58 class MappingP1
59 {
60 public:
61 
62  typedef arma::vec::fixed<dim+1> BaryPoint;
63  typedef arma::vec::fixed<spacedim> RealPoint;
64  typedef arma::mat::fixed<spacedim, dim+1> ElementMap;
65 
66  /**
67  * @brief Determines which additional quantities have to be computed.
68  *
69  * @param flags Update flags for required quantities.
70  * @return All necessary flags.
71  */
72  static UpdateFlags update_each(UpdateFlags flags);
73 
74  /**
75  * Map from reference element (barycentric coords) to global coord system.
76  * Matrix(3, dim+1) M: x_real = M * x_bary;
77  * M columns are real coordinates of nodes.
78  */
80 
81  /**
82  * Compute jacobian matrix for an element given by the @p coords element map.
83  */
84  static arma::mat::fixed<spacedim,dim> jacobian(const ElementMap &coords);
85 
86  /**
87  * Project given point in real coordinates to reference element (barycentic coordinates).
88  * Result vector have dimension dim()+1.
89  * Use RefElement<dim>::bary_to_local() to get local coordinates.
90  */
91  static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map);
92 
93  /**
94  * Project given point from reference element (barycentic coordinates) to real coordinates.
95  * Use RefElement<dim>::local_to_bary() to get barycentric coordinates in input.
96  */
97  static RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map);
98 
99  /**
100  * Clip a point given by barycentric cocordinates to the element.
101  * If the point is out of the element the closest point
102  * projection to the element surface is used.
103  */
104  static BaryPoint clip_to_element(BaryPoint &barycentric);
105 
106  /// Test if element contains given point.
107  static bool contains_point(arma::vec point, ElementAccessor<3> elm);
108 
109 
110 };
111 
112 
113 
114 
115 
116 
117 
118 
119 #endif /* MAPPING_P1_HH_ */
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:59
arma::vec::fixed< dim+1 > BaryPoint
Definition: mapping_p1.hh:62
arma::vec::fixed< spacedim > RealPoint
Definition: mapping_p1.hh:63
static BaryPoint clip_to_element(BaryPoint &barycentric)
Definition: mapping_p1.cc:91
static arma::mat::fixed< spacedim, dim > jacobian(const ElementMap &coords)
Definition: mapping_p1.cc:58
static bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
Definition: mapping_p1.cc:96
static RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map)
Definition: mapping_p1.cc:85
arma::mat::fixed< spacedim, dim+1 > ElementMap
Definition: mapping_p1.hh:64
static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map)
Definition: mapping_p1.cc:69
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:30
static const unsigned int dim_minus_one
Definition: mapping_p1.hh:38
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
ArmaVec< double, N > vec
Definition: armor.hh:933
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68