Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
nedelec_encode.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
37#ifndef __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_ENCODE__
38#define __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_ENCODE__
39
40#include "common/common.h"
42
44
48
49#ifdef UG_FOR_LUA
51#endif
52
53namespace ug{
54namespace Electromagnetism{
55
60template <typename TDomain>
62(
63 UserData<MathVector<TDomain::dim>, TDomain::dim> * function,
64 Edge * edge,
65 const typename TDomain::position_accessor_type& aaPos,
66 const int si,
67 number time
68)
69{
70 typedef typename TDomain::position_type position_type;
71
72 number dofValue;
73 MathVector<TDomain::dim> func_value;
74
75// Get the vector pointing along the edge from its beginning to its end:
76 position_type edgeVector = aaPos [edge->vertex(1)];
77 edgeVector -= aaPos [edge->vertex(0)];
78
79// Compute the quadrature
80
81 (*function) (func_value, aaPos [edge->vertex(0)], time, si);
82 dofValue = func_value * edgeVector;
83
84 (*function) (func_value, aaPos [edge->vertex(1)], time, si);
85 dofValue += func_value * edgeVector;
86
87 return dofValue / 2;
88}
89
94template <typename TGridFunction>
96(
97 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spFunction,
98 SmartPtr<TGridFunction> spGridFct,
99 size_t fct,
100 const SubsetGroup& ssGrp,
101 number time
102)
103{
104// Domain type, position_type and iterators
105 typedef typename TGridFunction::domain_type domain_type;
106 typedef typename TGridFunction::template traits<Edge>::const_iterator t_edge_iterator;
107
108// Get position accessor
109 const typename domain_type::position_accessor_type& aaPos = spGridFct->domain()->position_accessor();
110
111// Multiindex to access the components
112 std::vector<DoFIndex> ind;
113
114 UG_ASSERT (spFunction.valid (), "Invalid smart pointer to the user data.");
115
116// Loop over subsets
117 for (size_t i = 0; i < ssGrp.size (); i++)
118 {
119 // get subset index
120 const int si = ssGrp[i];
121
122 // skip if function is not defined in subset
123 if (! spGridFct->is_def_in_subset (fct, si))
124 {
125 UG_LOG ("ComputeNedelecDoFs warning: Function " << fct << " undefined in subdom " << si);
126 continue;
127 }
128
129 // Loop the edges
130 t_edge_iterator iter = spGridFct->template begin<Edge> (si);
131 t_edge_iterator iterEnd = spGridFct->template end<Edge> (si);
132 for (; iter != iterEnd; iter++)
133 {
134 Edge * pEdge = *iter;
135
136 // Get the multiindices
137 if (spGridFct->inner_dof_indices (pEdge, fct, ind) != 1)
138 UG_THROW ("ComputeNedelecDoFs:"
139 " More than one DoF per edge. Not the Nedelec-type-1 element?");
140
141 // Compute the value for the edge and write it into the grid function
142 DoFRef ((* spGridFct), ind[0])
143 = ComputeNedelecDoF<domain_type> (spFunction.get(), pEdge, aaPos, si, time);
144 }
145 }
146}
147
154template <typename TGridFunction>
156(
157 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spFunction,
158 SmartPtr<TGridFunction> spGridFct,
159 const char* cmp,
160 const char* subsets,
161 number time
162)
163{
164// Get function id by name
165 const size_t fct = spGridFct->fct_id_by_name (cmp);
166
167// Check that function exists
168 if (fct > spGridFct->num_fct())
169 UG_THROW("ComputeNedelecDoFs: Component with name '" << cmp << "' not found.");
170
171// Check the basis type
172 if (spGridFct->local_finite_element_id(fct).type () != LFEID::NEDELEC)
173 UG_THROW ("ComputeNedelecDoFs: Illegal grid function. It should be based on the Nedelec element.");
174
175// Create the subset group
176 SubsetGroup ssGrp (spGridFct->domain()->subset_handler ());
177 if (subsets != NULL)
178 ssGrp.add (TokenizeString (subsets));
179 else // add all subsets and remove lower dim subsets afterwards
180 ssGrp.add_all();
181
182// Compute the DoF values:
183 ComputeNedelecDoFs (spFunction, spGridFct, fct, ssGrp, time);
184}
185
186template <typename TGridFunction>
188(
189 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spFunction,
190 SmartPtr<TGridFunction> spGridFct,
191 const char* cmp,
192 number time
193)
194{
195 ComputeNedelecDoFs (spFunction, spGridFct, cmp, NULL, time);
196}
197
198template <typename TGridFunction>
200(
201 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spFunction,
202 SmartPtr<TGridFunction> spGridFct,
203 const char* cmp,
204 const char* subsets
205)
206{
207 ComputeNedelecDoFs (spFunction, spGridFct, cmp, subsets, 0.0);
208}
209
210template <typename TGridFunction>
212(
213 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spFunction,
214 SmartPtr<TGridFunction> spGridFct,
215 const char* cmp
216)
217{
218 ComputeNedelecDoFs (spFunction, spGridFct, cmp, NULL, 0.0);
219}
220
222
223#ifdef UG_FOR_LUA
224
230template <typename TGridFunction>
232(
233 const char* LuaFunction,
234 SmartPtr<TGridFunction> spGridFct,
235 const char* cmp,
236 const char* subsets,
237 number time
238)
239{
240 static const int dim = TGridFunction::dim;
241 SmartPtr<UserData<MathVector<dim>, dim> > spFunction
243 ComputeNedelecDoFs (spFunction, spGridFct, cmp, subsets, time);
244}
245template <typename TGridFunction>
247(
248 const char* LuaFunction,
249 SmartPtr<TGridFunction> spGridFct,
250 const char* cmp,
251 number time
252)
253{
254 ComputeNedelecDoFs (LuaFunction, spGridFct, cmp, NULL, time);
255}
256template <typename TGridFunction>
258(
259 const char* LuaFunction,
260 SmartPtr<TGridFunction> spGridFct,
261 const char* cmp,
262 const char* subsets
263)
264{
265 ComputeNedelecDoFs (LuaFunction, spGridFct, cmp, subsets, 0);
266}
267template <typename TGridFunction>
269(
270 const char* LuaFunction,
271 SmartPtr<TGridFunction> spGridFct,
272 const char* cmp
273)
274{
275 ComputeNedelecDoFs (LuaFunction, spGridFct, cmp, NULL, 0);
276}
277
278#endif // UG_FOR_LUA
279
284template <typename TPotGridFunc, typename TGridFunc>
286(
287 SmartPtr<TPotGridFunc> spPotGF,
289 size_t fct
290)
291{
292 typedef typename TGridFunc::template traits<Edge>::const_iterator t_edge_iterator;
293
294// Check the basis type
295 if (spGF->local_finite_element_id(fct).type () != LFEID::NEDELEC)
296 UG_THROW ("NedelecGradPotential: Illegal grid function. It should be based on the Nedelec element.");
297
298// Arrays for the indices in the vectors:
299 std::vector<size_t> vVertInd (1);
300 std::vector<DoFIndex> vEdgeInd (1);
301
302// Loop the edges
303 t_edge_iterator iterEnd = spGF->template end<Edge> ();
304 for (t_edge_iterator iter = spGF->template begin<Edge> (); iter != iterEnd; iter++)
305 {
306 number corner_val [2];
307 Edge * pEdge = *iter;
308
309 // Get the potential the ends of the edge:
310 for (size_t i = 0; i < 2; i++)
311 {
312 // Get the multiindices
313 if (spPotGF->inner_algebra_indices (pEdge->vertex(i), vVertInd) != 1)
314 UG_THROW ("NedelecGradPotential:"
315 " More than one DoF per vertex for the potential.");
316
317 // Set the Dirichlet entry
318 corner_val [i] = (* spPotGF) [vVertInd[0]];
319 }
320
321 // Compute the gradient:
322 if (spGF->inner_dof_indices (pEdge, fct, vEdgeInd) != 1)
323 UG_THROW ("NedelecGradPotential:"
324 " More than one DoF per edge. Not the Nedelec-Type-1 element?");
325 DoFRef ((* spGF), vEdgeInd[0])
326 -= corner_val [1] - corner_val [0];
327 }
328}
329
334template <typename TPotGridFunc, typename TGridFunc>
336(
337 SmartPtr<TPotGridFunc> spPotGF,
339 const char* cmp
340)
341{
342// Get function id by name
343 const size_t fct = spGF->fct_id_by_name (cmp);
344
345// Check that function exists
346 if (fct > spGF->num_fct())
347 UG_THROW("NedelecGradPotential: Component with name '" << cmp << "' not found.");
348
349// Distribute the potential
350 NedelecGradPotential (spPotGF, spGF, fct);
351}
352
357template <typename TPotGridFunc>
359(
360 SmartPtr<TPotGridFunc> spPotGF,
361 const SubsetGroup& ssGrp,
362 number value
363)
364{
365 typedef typename TPotGridFunc::template traits<Vertex>::const_iterator t_vert_iterator;
366
367// Array for the indices in the grid function:
368 std::vector<size_t> vVertInd (1);
369
370// Set the function to 0
371 spPotGF->set (0.0);
372
373// Loop over subsets
374 for (size_t i = 0; i < ssGrp.size (); i++)
375 {
376 // get subset index
377 const int si = ssGrp[i];
378
379 // Loop the vertices
380 t_vert_iterator iter = spPotGF->template begin<Vertex> (si);
381 t_vert_iterator iterEnd = spPotGF->template end<Vertex> (si);
382 for (; iter != iterEnd; iter++)
383 {
384 Vertex * pVertex = *iter;
385
386 // Get the multiindices
387 if (spPotGF->inner_algebra_indices (pVertex, vVertInd) != 1)
388 UG_THROW ("SubsetPotential:"
389 " More than one DoF per vertex for the potential.");
390
391 // Set the value
392 (* spPotGF) [vVertInd[0]] = value;
393 }
394 }
395}
396
401template <typename TPotGridFunc>
403(
404 SmartPtr<TPotGridFunc> spPotGF,
405 const char* subsets,
406 number value
407)
408{
409// Create the subset group
410 SubsetGroup ssGrp (spPotGF->domain()->subset_handler ());
411 if (subsets != NULL)
412 ssGrp.add (TokenizeString (subsets));
413 else // add all subsets and remove lower dim subsets afterwards
414 ssGrp.add_all();
415
416// Set the potential
417 SetSubsetVertVal (spPotGF, ssGrp, value);
418}
419
420} // end namespace Electromagnetism
421} // end namespace ug
422
423#endif // __H__UG__PLUGINS__ELECTROMAGNETISM__NEDELEC_ENCODE__
424
425/* End of File */
virtual Vertex * vertex(size_t index) const
void add(const char *name)
size_t size() const
vector< string > TokenizeString(const char *str, const char delimiter=',')
#define UG_ASSERT(expr, msg)
#define UG_THROW(msg)
#define UG_LOG(msg)
double number
void SetSubsetVertVal(SmartPtr< TPotGridFunc > spPotGF, const SubsetGroup &ssGrp, number value)
Definition nedelec_encode.h:359
void ComputeNedelecDoFs(SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, size_t fct, const SubsetGroup &ssGrp, number time)
Definition nedelec_encode.h:96
number ComputeNedelecDoF(UserData< MathVector< TDomain::dim >, TDomain::dim > *function, Edge *edge, const typename TDomain::position_accessor_type &aaPos, const int si, number time)
Definition nedelec_encode.h:62
void NedelecGradPotential(SmartPtr< TPotGridFunc > spPotGF, SmartPtr< TGridFunc > spGF, size_t fct)
Definition nedelec_encode.h:286
const number & DoFRef(const TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)