Flow123d  jenkins-Flow123d-linux-release-multijob-282
mapping_p1.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Class MappingP1 implements the affine transformation of
27  * the unit cell onto the actual cell.
28  * @author Jan Stebel
29  */
30 
31 #ifndef MAPPING_P1_HH_
32 #define MAPPING_P1_HH_
33 
34 #include "fem/mapping.hh"
35 
36 
37 /**
38  * Auxiliary templates to resolve nonzero matrix sizes for small dimensions.
39  * ( some compilers do not accept ternary operator in template parameters.)
40  */
41 template<unsigned int dim>
42 class MatrixSizes {
43 public:
44  static const unsigned int dim_minus_one = dim-1;
45 };
46 
47 template<>
48 class MatrixSizes<0> {
49 public:
50  static const unsigned int dim_minus_one = 0;
51 };
52 
53 
54 /**
55  * @brief Affine mapping between reference and actual cell.
56  *
57  * Class MappingP1 implements the affine transformation of
58  * the reference cell onto the actual cell.
59  *
60  * @param dim Dimension of the cells.
61  * @param spacedim Dimension of the Euclidean space.
62  */
63 template<unsigned int dim, unsigned int spacedim>
64 class MappingP1 : public Mapping<dim,spacedim>
65 {
66 public:
67 
68  /**
69  * @brief Constructor.
70  */
71  MappingP1();
72 
73  /**
74  * @brief Initializes the structures and computes static data.
75  *
76  * @param q Quadrature rule.
77  * @param flags Update flags.
78  * @return The computed mapping data.
79  */
81 
82  /**
83  * @brief Determines which additional quantities have to be computed.
84  *
85  * @param flags Update flags for required quantities.
86  * @return All necessary flags.
87  */
89 
90  /**
91  * @brief Calculates the mapping data on the actual cell.
92  *
93  * @param cell The actual cell.
94  * @param q Quadrature rule.
95  * @param data Precomputed mapping data.
96  * @param fv_data Data to be computed.
97  */
98  void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell,
99  const Quadrature<dim> &q,
100  MappingInternalData &data,
101  FEValuesData<dim,spacedim> &fv_data);
102 
103  /**
104  * @brief Calculates the mapping data on a side of a cell.
105  *
106  * @param cell The actual cell.
107  * @param sid Number of the side.
108  * @param q The quadrature rule with points on the side.
109  * @param data Precomputed mapping data.
110  * @param fv_data Data to be computed.
111  */
112  void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &cell,
113  unsigned int sid,
114  const Quadrature<dim> &q,
115  MappingInternalData &data,
116  FEValuesData<dim,spacedim> &fv_data);
117 
118 
119 private:
120 
121  /**
122  * @brief Auxiliary matrix of gradients of shape functions (used for
123  * computation of the Jacobian).
124  */
125  arma::mat::fixed<dim+1,dim> grad;
126 
127 };
128 
129 
130 
131 
132 
133 
134 
135 
136 #endif /* MAPPING_P1_HH_ */
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:78
MappingInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)
Initializes the structures and computes static data.
Definition: mapping_p1.cc:50
UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:95
static const unsigned int dim_minus_one
Definition: mapping_p1.hh:44
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:42
void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Calculates the mapping data on the actual cell.
Definition: mapping_p1.cc:112
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:44
void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &cell, unsigned int sid, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Calculates the mapping data on a side of a cell.
Definition: mapping_p1.cc:188
MappingP1()
Constructor.
Definition: mapping_p1.cc:45
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:64
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:125
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:57
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:110