Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
nedelec_gf_user_data.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2014: 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
36#ifndef __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_GF_USER_DATA__
37#define __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_GF_USER_DATA__
38
39#include <map>
40#include <vector>
41
42// ug4 headers
43#include "common/common.h"
44#include "lib_grid/lg_base.h"
48
49// Nedelec-type-1 headers:
50#include "nedelec_local_ass.h"
51
52#include "em_material.h"
53
54namespace ug{
55namespace Electromagnetism{
56
61template <typename TGridFunc>
64 <NedelecGridFunctionData<TGridFunc>, MathVector<TGridFunc::dim>, TGridFunc::dim>
65{
66public:
68 typedef typename TGridFunc::domain_type domain_type;
69
71 static const int dim = domain_type::dim;
72
74 typedef typename domain_type::position_type position_type;
75
76private:
79
81 size_t m_fct;
82
83public:
86 (
87 SmartPtr<TGridFunc> spGridFct,
88 const char * cmp
89 )
90 : m_spGF (spGridFct)
91 {
92 // Get function id by name:
93 try
94 {
95 m_fct = m_spGF->fct_id_by_name (cmp);
96 }
97 UG_CATCH_THROW ("NedelecGridFunctionData: Function space does not contain"
98 " a function with name " << cmp << ".");
99
100 // Get local finite element id and check whether this is the Nedelec element:
101 if (m_spGF->local_finite_element_id(m_fct).type () != LFEID::NEDELEC)
102 UG_THROW ("NedelecGridFunctionData: Function " << cmp
103 << "is not based on the Nedelec element.");
104 };
105
107 virtual bool continuous () const {return false;}
108
110 virtual bool requires_grid_fct () const {return true;}
111
113 template <int refDim>
115 (
116 MathVector<dim> vValue[],
117 const MathVector<dim> vGlobIP[],
118 number time,
119 int si,
120 GridObject * elem,
121 const MathVector<dim> vCornerCoords[],
122 const MathVector<refDim> vLocIP[],
123 const size_t nip,
124 LocalVector * u,
125 bool bDeriv,
126 int s,
127 std::vector<std::vector<MathVector<dim> > > vvvDeriv[],
128 const MathMatrix<refDim, dim> * vJT = NULL
129 ) const
130 {
131 static const size_t maxEdges = NedelecInterpolation<domain_type, refDim>::maxNumEdges;
132
133 // Derivatives are not implemented
134 if (bDeriv)
135 UG_THROW ("NedelecGridFunctionData: Derivatives are not implemented.");
136
137 // Get multiindices of element
138 std::vector<DoFIndex> ind;
139 m_spGF->dof_indices (elem, m_fct, ind);
140 if (ind.size() > maxEdges)
141 UG_THROW ("NedelecGridFunctionData: Illegal number of DoFs!");
142
143 // The DoF values of the grid function
144 number dofValues [maxEdges];
145 for (size_t sh = 0; sh < ind.size(); ++sh)
146 dofValues[sh] = DoFRef (*m_spGF, ind[sh]);
147
148 // Compute the values
150 (m_spGF->domain().get(), elem, vCornerCoords, dofValues, vLocIP, nip, vValue);
151 };
152};
153
158template <typename TGridFunc>
160 : public StdDependentUserData
161 <NedelecCurlData<TGridFunc>, MathVector<TGridFunc::dim>, TGridFunc::dim>
162{
163public:
165 typedef typename TGridFunc::domain_type domain_type;
166
168 static const int dim = domain_type::dim;
169
171 typedef typename domain_type::position_type position_type;
172
173private:
176
178 size_t m_fct;
179
180public:
183 (
184 SmartPtr<TGridFunc> spGridFct,
185 const char * cmp
186 )
187 : m_spGF (spGridFct)
188 {
189 // Get function id by name:
190 try
191 {
192 m_fct = m_spGF->fct_id_by_name (cmp);
193 }
194 UG_CATCH_THROW ("NedelecCurlData: Function space does not contain"
195 " a function with name " << cmp << ".");
196
197 // Get local finite element id and check whether this is the Nedelec element:
198 if (m_spGF->local_finite_element_id(m_fct).type () != LFEID::NEDELEC)
199 UG_THROW ("NedelecCurlData: Function " << cmp
200 << "is not based on the Nedelec element.");
201 };
202
204 virtual bool continuous () const {return false;}
205
207 virtual bool requires_grid_fct () const {return true;}
208
210 template <int refDim>
212 (
213 MathVector<dim> vValue[],
214 const MathVector<dim> vGlobIP[],
215 number time,
216 int si,
217 GridObject * elem,
218 const MathVector<dim> vCornerCoords[],
219 const MathVector<refDim> vLocIP[],
220 const size_t nip,
221 LocalVector * u,
222 bool bDeriv,
223 int s,
224 std::vector<std::vector<MathVector<dim> > > vvvDeriv[],
225 const MathMatrix<refDim, dim> * vJT = NULL
226 ) const
227 {
228 static const size_t maxEdges = NedelecInterpolation<domain_type, refDim>::maxNumEdges;
229
230 // Derivatives are not implemented
231 if (bDeriv)
232 UG_THROW ("NedelecCurlData: Derivatives are not implemented.");
233
234 // Get multiindices of element
235 std::vector<DoFIndex> ind;
236 m_spGF->dof_indices (elem, m_fct, ind);
237 if (ind.size() > maxEdges)
238 UG_THROW ("NedelecGridFunctionData: Illegal number of DoFs!");
239
240 // The DoF values of the grid function
241 number dofValues [maxEdges];
242 for (size_t sh = 0; sh < ind.size(); ++sh)
243 dofValues[sh] = DoFRef (*m_spGF, ind[sh]);
244
245 // Compute the curl
246 MathVector<dim> curl;
248 (m_spGF->domain().get(), elem, vCornerCoords, dofValues, curl);
249 for (size_t ip = 0; ip < nip; ip++)
250 vValue[ip] = curl;
251 };
252};
253
261template <typename TGridFunc>
263 : public StdDependentUserData
264 <NedelecSigmaEData<TGridFunc>, MathVector<TGridFunc::dim>, TGridFunc::dim>
265{
266public:
268 typedef typename TGridFunc::domain_type domain_type;
269
271 static const int dim = domain_type::dim;
272
274 typedef typename domain_type::position_type position_type;
275
276private:
279
281 size_t m_fct;
282
285
286public:
289 (
290 SmartPtr<TGridFunc> spGridFct,
291 const char * cmp,
292 SmartPtr<EMaterial<domain_type> > emMatherial
293 )
294 : m_spGF (spGridFct), m_spEMaterial (emMatherial)
295 {
296 // Get function id by name:
297 try
298 {
299 m_fct = m_spGF->fct_id_by_name (cmp);
300 }
301 UG_CATCH_THROW ("NedelecGridFunctionData: Function space does not contain"
302 " a function with name " << cmp << ".");
303
304 // Get local finite element id and check whether this is the Nedelec element:
305 if (m_spGF->local_finite_element_id(m_fct).type () != LFEID::NEDELEC)
306 UG_THROW ("NedelecGridFunctionData: Function " << cmp
307 << "is not based on the Nedelec element.");
308 };
309
311 virtual bool continuous () const {return false;}
312
314 virtual bool requires_grid_fct () const {return true;}
315
317 template <int refDim>
319 (
320 MathVector<dim> vValue[],
321 const MathVector<dim> vGlobIP[],
322 number time,
323 int si,
324 GridObject * elem,
325 const MathVector<dim> vCornerCoords[],
326 const MathVector<refDim> vLocIP[],
327 const size_t nip,
328 LocalVector * u,
329 bool bDeriv,
330 int s,
331 std::vector<std::vector<MathVector<dim> > > vvvDeriv[],
332 const MathMatrix<refDim, dim> * vJT = NULL
333 ) const
334 {
335 static const size_t maxEdges = NedelecInterpolation<domain_type, refDim>::maxNumEdges;
336
337 // Derivatives are not implemented
338 if (bDeriv)
339 UG_THROW ("NedelecSigmaEData: Derivatives are not implemented.");
340
341 // Get the material data
342 number mu, sigma;
343 if (m_spEMaterial->get_mu_sigma (si, mu, sigma))
344 UG_THROW ("NedelecSigmaEData: No material data set to subset" << si
345 << " (or this is a low-dim. domain).");
346 if (sigma == 0) // no currents in insulators
347 {
348 for (size_t i = 0; i < nip; i++) vValue [i] = 0.0;
349 return;
350 }
351
352 // Get multiindices of element
353 std::vector<DoFIndex> ind;
354 m_spGF->dof_indices (elem, m_fct, ind);
355 if (ind.size() > maxEdges)
356 UG_THROW ("NedelecGridFunctionData: Illegal number of DoFs!");
357
358 // The DoF values of the grid function
359 number dofValues [maxEdges];
360 for (size_t sh = 0; sh < ind.size(); ++sh)
361 dofValues[sh] = DoFRef (*m_spGF, ind[sh]);
362
363 // Compute the values of the grid function
365 (m_spGF->domain().get(), elem, vCornerCoords, dofValues, vLocIP, nip, vValue);
366
367 // Multiply the values by \f$ \sigma \f$
368 for (size_t i = 0; i < nip; i++)
369 vValue [i] *= sigma;
370 };
371};
372
373} // end namespace Electromagnetism
374} // end namespace ug
375
376#endif // __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_GF_USER_DATA__
377
378/* End of File */
parameterString s
Definition Biogas.lua:2
T * get()
Class for subdomain-dependent data for the E-based formulated problems.
Definition em_material.h:63
Definition nedelec_gf_user_data.h:162
TGridFunc::domain_type domain_type
Type of domain.
Definition nedelec_gf_user_data.h:165
size_t m_fct
component of function
Definition nedelec_gf_user_data.h:178
domain_type::position_type position_type
Type of position coordinates (e.g. position_type)
Definition nedelec_gf_user_data.h:171
virtual bool requires_grid_fct() const
Returns true to get the grid element in the evaluation routine.
Definition nedelec_gf_user_data.h:207
static const int dim
World dimension.
Definition nedelec_gf_user_data.h:168
SmartPtr< TGridFunc > m_spGF
grid function
Definition nedelec_gf_user_data.h:175
NedelecCurlData(SmartPtr< TGridFunc > spGridFct, const char *cmp)
constructor
Definition nedelec_gf_user_data.h:183
void eval_and_deriv(MathVector< dim > vValue[], const MathVector< dim > vGlobIP[], number time, int si, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< refDim > vLocIP[], const size_t nip, LocalVector *u, bool bDeriv, int s, std::vector< std::vector< MathVector< dim > > > vvvDeriv[], const MathMatrix< refDim, dim > *vJT=NULL) const
Performs the main computations:
Definition nedelec_gf_user_data.h:212
virtual bool continuous() const
The vector field retrieved from the Nedelec-type 1 (Whitney-1) dofs are not continuous.
Definition nedelec_gf_user_data.h:204
Definition nedelec_gf_user_data.h:65
static const int dim
World dimension.
Definition nedelec_gf_user_data.h:71
TGridFunc::domain_type domain_type
Type of domain.
Definition nedelec_gf_user_data.h:68
virtual bool continuous() const
The vector field retrieved from the Nedelec-type 1 (Whitney-1) dofs are not continuous.
Definition nedelec_gf_user_data.h:107
virtual bool requires_grid_fct() const
Returns true to get the grid element in the evaluation routine.
Definition nedelec_gf_user_data.h:110
SmartPtr< TGridFunc > m_spGF
grid function
Definition nedelec_gf_user_data.h:78
domain_type::position_type position_type
Type of position coordinates (e.g. position_type)
Definition nedelec_gf_user_data.h:74
NedelecGridFunctionData(SmartPtr< TGridFunc > spGridFct, const char *cmp)
constructor
Definition nedelec_gf_user_data.h:86
size_t m_fct
component of function
Definition nedelec_gf_user_data.h:81
void eval_and_deriv(MathVector< dim > vValue[], const MathVector< dim > vGlobIP[], number time, int si, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< refDim > vLocIP[], const size_t nip, LocalVector *u, bool bDeriv, int s, std::vector< std::vector< MathVector< dim > > > vvvDeriv[], const MathMatrix< refDim, dim > *vJT=NULL) const
Performs the main computations:
Definition nedelec_gf_user_data.h:115
Interpolation of the Nedelec dofs and their curls.
Definition nedelec_local_ass.h:399
static void curl(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], MathVector< WDim > &curl_vec)
computes curl of the function
Definition nedelec_local_ass.h:434
static void value(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< refDim > local[], const size_t n_pnt, MathVector< WDim > values[])
computes the values at given points
Definition nedelec_local_ass.h:411
Definition nedelec_gf_user_data.h:265
TGridFunc::domain_type domain_type
Type of domain.
Definition nedelec_gf_user_data.h:268
SmartPtr< TGridFunc > m_spGF
grid function
Definition nedelec_gf_user_data.h:278
NedelecSigmaEData(SmartPtr< TGridFunc > spGridFct, const char *cmp, SmartPtr< EMaterial< domain_type > > emMatherial)
constructor
Definition nedelec_gf_user_data.h:289
static const int dim
World dimension.
Definition nedelec_gf_user_data.h:271
SmartPtr< EMaterial< domain_type > > m_spEMaterial
subdomain-dependent data (propertiels of the materials)
Definition nedelec_gf_user_data.h:284
size_t m_fct
component of function
Definition nedelec_gf_user_data.h:281
virtual bool requires_grid_fct() const
Returns true to get the grid element in the evaluation routine.
Definition nedelec_gf_user_data.h:314
virtual bool continuous() const
The vector field retrieved from the Nedelec-type 1 (Whitney-1) dofs are not continuous.
Definition nedelec_gf_user_data.h:311
domain_type::position_type position_type
Type of position coordinates (e.g. position_type)
Definition nedelec_gf_user_data.h:274
void eval_and_deriv(MathVector< dim > vValue[], const MathVector< dim > vGlobIP[], number time, int si, GridObject *elem, const MathVector< dim > vCornerCoords[], const MathVector< refDim > vLocIP[], const size_t nip, LocalVector *u, bool bDeriv, int s, std::vector< std::vector< MathVector< dim > > > vvvDeriv[], const MathMatrix< refDim, dim > *vJT=NULL) const
Performs the main computations:
Definition nedelec_gf_user_data.h:319
number time() const
const MathVector< dim > & ip(size_t s, size_t ip) const
#define UG_CATCH_THROW(msg)
#define UG_THROW(msg)
double number
const number & DoFRef(const TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)