Line data Source code
1 : /*
2 : * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Andreas Vogel
4 : *
5 : * This file is part of UG4.
6 : *
7 : * UG4 is free software: you can redistribute it and/or modify it under the
8 : * terms of the GNU Lesser General Public License version 3 (as published by the
9 : * Free Software Foundation) with the following additional attribution
10 : * requirements (according to LGPL/GPL v3 §7):
11 : *
12 : * (1) The following notice must be displayed in the Appropriate Legal Notices
13 : * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 : *
15 : * (2) The following notice must be displayed at a prominent place in the
16 : * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 : *
18 : * (3) The following bibliography is recommended for citation and must be
19 : * preserved in all covered files:
20 : * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 : * parallel geometric multigrid solver on hierarchically distributed grids.
22 : * Computing and visualization in science 16, 4 (2013), 151-164"
23 : * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 : * flexible software system for simulating pde based models on high performance
25 : * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 : *
27 : * This program is distributed in the hope that it will be useful,
28 : * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 : * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 : * GNU Lesser General Public License for more details.
31 : */
32 :
33 : #ifndef __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__LOCAL_FINITE_ELEMENT_ID__
34 : #define __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__LOCAL_FINITE_ELEMENT_ID__
35 :
36 : #include <sstream>
37 :
38 : namespace ug{
39 :
40 :
41 : // Doxygen group
42 : ///////////////////////////////////////////////////////////////////////////////
43 : /**
44 : * \brief provides Local Finite Elements.
45 : *
46 : * The Local Finite Element section is used to describe finite element spaces
47 : * by their definition on reference elements.
48 : *
49 : * A Finite Element is defined as a triplet \f$ \{ K, P, \Sigma \} \f$
50 : * (See e.g.
51 : * Ciarlet, P., "Basis Error Estimates for Elliptic Problems", North-Holland,
52 : * Amsterdam, 1991, p. 93;
53 : * or Ern, A. and Guermond J.L., "Theory and Practice of Finite Elements",
54 : * Springer, 2004, p. 19), where
55 : *
56 : * <ol>
57 : * <li> \f$ K \f$ is a compact, connected, Lipschitz subset of \f$\mathbb{R}^d\f$
58 : * with non-empty interior
59 : * <li> \f$ P \f$ is a vector space of functions \f$p: K \mapsto \mathbb{R}^m \f$
60 : * with an integer \f$m > 0 \f$ (usually \f$m=1\f$ or \f$m=d\f$)
61 : * <li> \f$\Sigma\f$ is a set of \f$ n_{sh} \f$ linear forms \f$ \sigma_1, \dots,
62 : * \sigma_{n_{sh}} \f$ acting on the elements of \f$ P \f$, such that the
63 : * linear mapping
64 : * \f[
65 : * p \mapsto ( \sigma_1(p), \dots, \sigma_{n_{sh}}(p)) \in \mathbb{R}^{n_{sh}}
66 : * \f]
67 : * is bijective. These linear forms are called local degrees of freedom.
68 : * </ol>
69 : *
70 : * Since the mapping is bijective, there exist a basis \f$\{\phi_1, \dots,
71 : * \phi_{n_{sh}} \subset P\f$ such that
72 : * \f[
73 : * \sigma_i (\phi_j) = \delta_{ij}, \qquad 1 \leq i,j \leq n_{sh}.
74 : * \f]
75 : *
76 : * This set is called the set of local shape functions. The implemented
77 : * counterpart is the class LocalShapeFunctionSet.
78 : *
79 : * The set of local degrees of freedom finds its counterpart in the class
80 : * ILocalDoFSet.
81 : *
82 : * \defgroup lib_disc_local_finite_elements Local Finite Elements
83 : * \ingroup lib_discretization
84 : */
85 : ///////////////////////////////////////////////////////////////////////////////
86 :
87 : /// \ingroup lib_disc_local_finite_elements
88 : /// @{
89 :
90 : /// Identifier for Local Finite Elements
91 : /**
92 : * This Class is used to distinguish between different local finite elements (lfe).
93 : * Each lfe has a unique Space Type (e.g. Lagrange, DG) and the order of trial
94 : * functions. It the function space is p-adaptive, the enum ADAPTIVE is set as
95 : * order.
96 : */
97 : class LFEID
98 : {
99 : public:
100 : /// Space Type
101 : enum SpaceType
102 : {
103 : NONE = -1,
104 : LAGRANGE = 0,
105 : CROUZEIX_RAVIART,
106 : PIECEWISE_CONSTANT,
107 : DG,
108 : MINI,
109 : NEDELEC,
110 : USER_DEFINED,
111 : NUM_SPACE_TYPES
112 : };
113 :
114 : /// special possibilities for order
115 : enum {ADAPTIV = -1, INVALID = -10};
116 :
117 : public:
118 : /// default constructor
119 0 : LFEID() : m_type(NONE), m_dim(-1), m_order(INVALID) {}
120 :
121 : /// constructor with values
122 : LFEID(SpaceType type, int dim, int order)
123 0 : : m_type(type), m_dim(dim), m_order(order) {}
124 :
125 : /// returns the order of the local finite element
126 0 : int order() const {return m_order;}
127 :
128 : /// returns the space dimension of the local finite element
129 0 : int dim() const {return m_dim;}
130 :
131 : /// returns the type of the local finite element
132 0 : SpaceType type() const {return m_type;}
133 :
134 : /// equality check
135 : bool operator==(const LFEID& v) const
136 : {
137 0 : return (m_type == v.m_type && m_dim==v.m_dim && m_order==v.m_order);
138 : }
139 :
140 : /// inequality check
141 : bool operator!=(const LFEID& v) const {return !((*this)==v);}
142 :
143 : /// operator <
144 : /** returns comparison which set id is regarded lesser
145 : * The Local Finite Elements are ordered by type first, then by dimension and
146 : * then by increasing order.
147 : */
148 : bool operator<(const LFEID& v) const
149 : {
150 0 : if(m_type != v.m_type) return m_type < v.m_type;
151 0 : else if(m_dim != v.m_dim) return m_dim < v.m_dim;
152 0 : else return m_order < v.m_order;
153 : }
154 :
155 : /// operator >
156 : bool operator>(const LFEID& v) const
157 : {
158 : if(m_type != v.m_type) return m_type > v.m_type;
159 : else if(m_dim != v.m_dim) return m_dim > v.m_dim;
160 : else return m_order > v.m_order;
161 : }
162 :
163 : /// operator <=
164 : bool operator<=(const LFEID& v) const
165 : {
166 : return (*this < v || *this == v);
167 : }
168 :
169 : /// operator >=
170 : bool operator>=(const LFEID& v) const
171 : {
172 : return (*this > v || *this == v);
173 : }
174 :
175 : friend std::ostream& operator<<(std::ostream& out, const LFEID& v);
176 :
177 : private:
178 : /// Space type
179 : SpaceType m_type;
180 :
181 : /// dimension
182 : int m_dim;
183 :
184 : /// Order
185 : int m_order;
186 : };
187 :
188 : /// writes the Identifier to the output stream
189 : std::ostream& operator<<(std::ostream& out, const LFEID& v);
190 :
191 : /// returns the LFEID for a combination of Space and order
192 : LFEID ConvertStringToLFEID(const char* type, int dim, int order);
193 :
194 : /// returns the LFEID
195 : LFEID ConvertStringToLFEID(const char* type, int dim);
196 :
197 : /// @}
198 :
199 : } // end namespace ug
200 :
201 : #endif /* __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__LOCAL_FINITE_ELEMENT_ID__ */
|