Flow123d  release_2.2.0-37-g336ee74
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 
25 
26 /**
27  * Auxiliary templates to resolve nonzero matrix sizes for small dimensions.
28  * ( some compilers do not accept ternary operator in template parameters.)
29  */
30 template<unsigned int dim>
31 class MatrixSizes {
32 public:
33  static const unsigned int dim_minus_one = dim-1;
34 };
35 
36 template<>
37 class MatrixSizes<0> {
38 public:
39  static const unsigned int dim_minus_one = 0;
40 };
41 
42 
43 /**
44  * @brief Affine mapping between reference and actual cell.
45  *
46  * Class MappingP1 implements the affine transformation of
47  * the reference cell onto the actual cell.
48  *
49  * @param dim Dimension of the cells.
50  * @param spacedim Dimension of the Euclidean space.
51  */
52 template<unsigned int dim, unsigned int spacedim>
53 class MappingP1 : public Mapping<dim,spacedim>
54 {
55 public:
56 
57  /**
58  * @brief Constructor.
59  */
60  MappingP1();
61 
62  /**
63  * @brief Initializes the structures and computes static data.
64  *
65  * @param q Quadrature rule.
66  * @param flags Update flags.
67  * @return The computed mapping data.
68  */
69  MappingInternalData *initialize(const Quadrature<dim> &q, UpdateFlags flags);
70 
71  /**
72  * @brief Determines which additional quantities have to be computed.
73  *
74  * @param flags Update flags for required quantities.
75  * @return All necessary flags.
76  */
77  UpdateFlags update_each(UpdateFlags flags);
78 
79  /**
80  * @brief Calculates the mapping data on the actual cell.
81  *
82  * @param cell The actual cell.
83  * @param q Quadrature rule.
84  * @param data Precomputed mapping data.
85  * @param fv_data Data to be computed.
86  */
87  void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell,
88  const Quadrature<dim> &q,
89  MappingInternalData &data,
91 
92  /**
93  * @brief Calculates the mapping data on a side of a cell.
94  *
95  * @param cell The actual cell.
96  * @param sid Number of the side.
97  * @param q The quadrature rule with points on the side.
98  * @param data Precomputed mapping data.
99  * @param fv_data Data to be computed.
100  */
101  void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &cell,
102  unsigned int sid,
103  const Quadrature<dim> &q,
104  MappingInternalData &data,
105  FEValuesData<dim,spacedim> &fv_data);
106 
107 
108 private:
109 
110  /**
111  * @brief Auxiliary matrix of gradients of shape functions (used for
112  * computation of the Jacobian).
113  */
114  arma::mat::fixed<dim+1,dim> grad;
115 
116 };
117 
118 
119 
120 
121 
122 
123 
124 
125 #endif /* MAPPING_P1_HH_ */
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:33
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:31
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:33
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:53
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:114
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:46
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:99