Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
nedelec.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3 * Author: Dmitry Logashenko
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/*
34 * This file contains implementations of the local shape function set for
35 * the so-called Nedelec (or Whitney-1) elements.
36 */
37#ifndef __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__NEDELEC__NEDELEC__
38#define __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__NEDELEC__NEDELEC__
39
41#include "../local_dof_set.h"
43
44namespace ug{
45
47// Nedelec Set
49
51template <typename TRefElem>
53{
54 protected:
56 static const int refDim = TRefElem::dim;
57
58 public:
61 {
62 if(refDim < 2)
63 {
64 // No dofs if the dimension is less than 2:
65 nsh = 0;
66 m_vLocalDoF.clear();
67 return;
68 }
69
70 const TRefElem& rRefElem = Provider<TRefElem>::get();
71
72 nsh = rRefElem.num(1); // number of the edges
73 // set local DoFs (all located at the edges)
74 m_vLocalDoF.resize(nsh);
75 for(size_t i = 0; i < nsh; ++i)
76 m_vLocalDoF[i] = LocalDoF(1, i, 0);
77 }
78
80 ReferenceObjectID roid() const {return TRefElem::REFERENCE_OBJECT_ID;}
81
83 size_t num_dof() const {return nsh;};
84
86 size_t num_dof(ReferenceObjectID type) const
87 {
88 if(ReferenceElementDimension(type) == 1) return 1;
89 else return 0;
90 }
91
93 const LocalDoF& local_dof(size_t dof) const {return m_vLocalDoF[dof];}
94
96 bool exact_position_available() const {return true;};
97
98 protected:
100 size_t nsh;
101
103 std::vector<LocalDoF> m_vLocalDoF;
104};
105
112template <typename TRefElement>
114
116template <>
118: public NedelecLDS<ReferenceTriangle>,
119 public
121 <
122 NedelecLSFS<ReferenceTriangle>,
123 ReferenceTriangle::dim,
124 MathVector<ReferenceTriangle::dim>,
125 MathMatrix<ReferenceTriangle::dim, ReferenceTriangle::dim>
126 >
127{
128 public:
131
133 static const size_t order = 1;
134
136 static const int dim = reference_element_type::dim;
137
138 private:
141
142 public:
144 typedef base_type::shape_type shape_type;
145
147 typedef base_type::grad_type grad_type;
148
149 protected:
151 static const size_t nsh = reference_element_type::numEdges;
152
153 public:
156
158 inline bool continuous() const {return false;}
159
161 inline size_t num_sh() const {return nsh;}
162
164 inline bool position(size_t i, MathVector<dim>& pos) const
165 {
166 switch(i)
167 {
168 case 0: pos[0] = 0.5;
169 pos[1] = 0.0;
170 return true;
171 case 1: pos[0] = 0.5;
172 pos[1] = 0.5;
173 return true;
174 case 2: pos[0] = 0.0;
175 pos[1] = 0.5;
176 return true;
177 default: UG_THROW("NedelecLSFS: shape function "<<i<<
178 " not found. Only "<<nsh<<" shapes present.");
179 }
180 }
181
183 inline MathVector<dim> shape(const size_t i, const MathVector<dim>& x) const
184 {
185 UG_THROW ("NedelecLSFS: Nedelec shapes cannot be computed in the reference space.");
186 }
187
189 inline void grad(MathMatrix<dim,dim>& g, const size_t i, const MathVector<dim>& x) const
190 {
191 UG_THROW ("NedelecLSFS: Gradients of the Nedelec shapes cannot be computed in the reference space.");
192 }
193};
194
196template <>
198: public NedelecLDS<ReferenceTetrahedron>,
199 public
201 <
202 NedelecLSFS<ReferenceTetrahedron>,
203 ReferenceTetrahedron::dim,
204 MathVector<ReferenceTetrahedron::dim>,
205 MathMatrix<ReferenceTetrahedron::dim, ReferenceTetrahedron::dim>
206 >
207{
208 public:
211
213 static const size_t order = 1;
214
216 static const int dim = reference_element_type::dim;
217
218 private:
221
222 public:
224 typedef base_type::shape_type shape_type;
225
227 typedef base_type::grad_type grad_type;
228
229 protected:
231 static const size_t nsh = reference_element_type::numEdges;
232
233 public:
236
238 inline bool continuous() const {return false;}
239
241 inline size_t num_sh() const {return nsh;}
242
244 inline bool position(size_t i, MathVector<dim>& pos) const
245 {
246 switch(i)
247 {
248 case 0: pos[0] = 0.5;
249 pos[1] = 0.0;
250 pos[2] = 0.0;
251 return true;
252 case 1: pos[0] = 0.5;
253 pos[1] = 0.5;
254 pos[2] = 0.0;
255 return true;
256 case 2: pos[0] = 0.0;
257 pos[1] = 0.5;
258 pos[2] = 0.0;
259 return true;
260 case 3: pos[0] = 0.0;
261 pos[1] = 0.0;
262 pos[2] = 0.5;
263 return true;
264 case 4: pos[0] = 0.5;
265 pos[1] = 0.0;
266 pos[2] = 0.5;
267 return true;
268 case 5: pos[0] = 0.0;
269 pos[1] = 0.5;
270 pos[2] = 0.5;
271 return true;
272 default: UG_THROW("NedelecLSFS: shape function "<<i<<
273 " not found. Only "<<nsh<<" shapes present.");
274 }
275 }
276
278 inline MathVector<dim> shape(const size_t i, const MathVector<dim>& x) const
279 {
280 UG_THROW ("NedelecLSFS: Nedelec shapes cannot be computed in the reference space.");
281 }
282
284 inline void grad(MathMatrix<dim,dim>& g, const size_t i, const MathVector<dim>& x) const
285 {
286 UG_THROW ("NedelecLSFS: Gradients of the Nedelec shapes cannot be computed in the reference space.");
287 }
288};
289
290} // namespace ug
291
292#endif // __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__NEDELEC__NEDELEC__
293
294/* End of File */
static interface for trial spaces
Definition local_shape_function_set.h:176
Definition local_dof_set.h:58
A class for fixed size, dense matrices.
Definition math_matrix.h:63
a mathematical Vector with N entries.
Definition math_vector.h:97
Nedelec, i.e. the edge local dof set.
Definition nedelec.h:53
size_t nsh
number of shapes (== number of edges)
Definition nedelec.h:100
size_t num_dof(ReferenceObjectID type) const
returns the number of DoFs on a sub-geometric object type
Definition nedelec.h:86
const LocalDoF & local_dof(size_t dof) const
returns the dof storage
Definition nedelec.h:93
bool exact_position_available() const
returns if the local dof position are exact
Definition nedelec.h:96
size_t num_dof() const
returns the total number of DoFs on the finite element
Definition nedelec.h:83
std::vector< LocalDoF > m_vLocalDoF
association to elements
Definition nedelec.h:103
static const int refDim
dimension of reference element
Definition nedelec.h:56
ReferenceObjectID roid() const
returns the type of reference element
Definition nedelec.h:80
NedelecLDS()
constructor
Definition nedelec.h:60
base_type::grad_type grad_type
Gradient type.
Definition nedelec.h:227
void grad(MathMatrix< dim, dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition nedelec.h:284
MathVector< dim > shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition nedelec.h:278
base_type::shape_type shape_type
Shape type.
Definition nedelec.h:224
NedelecLSFS()
Constructor.
Definition nedelec.h:235
ReferenceTetrahedron reference_element_type
Reference Element type.
Definition nedelec.h:210
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition nedelec.h:244
size_t num_sh() const
Definition nedelec.h:241
bool continuous() const
returns if space constructs continuous functions
Definition nedelec.h:238
BaseLSFS< NedelecLSFS< reference_element_type >, dim, MathVector< dim >, MathMatrix< dim, dim > > base_type
Base class.
Definition nedelec.h:220
NedelecLSFS()
Constructor.
Definition nedelec.h:155
MathVector< dim > shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition nedelec.h:183
size_t num_sh() const
Definition nedelec.h:161
ReferenceTriangle reference_element_type
Reference Element type.
Definition nedelec.h:130
base_type::grad_type grad_type
Gradient type.
Definition nedelec.h:147
bool continuous() const
returns if space constructs continuous functions
Definition nedelec.h:158
BaseLSFS< NedelecLSFS< reference_element_type >, dim, MathVector< dim >, MathMatrix< dim, dim > > base_type
Base class.
Definition nedelec.h:140
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition nedelec.h:164
base_type::shape_type shape_type
Shape type.
Definition nedelec.h:144
void grad(MathMatrix< dim, dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition nedelec.h:189
Definition nedelec.h:113
static TClass & get()
returns a singleton based on the identifier
Definition provider.h:54
Definition reference_element.h:494
Definition reference_element.h:394
#define UG_THROW(msg)
Definition error.h:57
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
int ReferenceElementDimension(ReferenceObjectID roid)
returns the reference element dimension at run-time
Definition reference_element_util.h:40