Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
nedelec_dirichlet_impl.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
33/*
34 * Implementation of the functions for the Dirichlet boundary conditions
35 * for the Nedelec-based discretizations of the rot-rot operator.
36 */
37
38#ifdef UG_FOR_LUA
40#endif
41
42namespace ug{
43namespace Electromagnetism{
44
48template <typename TDomain, typename TAlgebra>
50(
51 const char* f_name,
52 const char* ss_names
53)
54{
55 typedef std::vector<std::string> t_f_vec;
56
57 std::vector<std::string> subsets;
58 TokenizeString (std::string (ss_names), subsets);
59
60// Check the data
61 std::string s_f_name (f_name);
62 size_t n_func;
63 if ((n_func = TokenizeString(s_f_name).size ()) > 1)
64 UG_THROW("NedelecDirichletBC:"
65 " Only single functions function allowed in the specification of"
66 " Dirichlet values, but the following functions given:"
67 << f_name);
68
69// Loop the subsets
70 for (size_t i = 0; i < subsets.size (); i++)
71 {
72 // If the map entry does not exist, create and initialize with all the functions
73 if (m_mDirichletSS.find (subsets [i]) == m_mDirichletSS.end ())
74 m_mDirichletSS [subsets [i]] = m_vDirichletFunc;
75 // Exclude the specified function from the list of the implicit specifications
76 std::vector<std::string> & imp_funcs = m_mDirichletSS [subsets [i]];
77 if (n_func == 0) continue;
78 t_f_vec::iterator f_entry = std::find (imp_funcs.begin (), imp_funcs.end (), s_f_name);
79 if (f_entry == imp_funcs.end ())
80 UG_THROW ("NedelecDirichletBC: Wrong speficifation of the Dirichlet BC for subset "
81 << subsets [i] << ": func " << f_name <<
82 " not present in the object or is specified twice for the subsets.");
83 imp_funcs.erase (f_entry);
84 }
85}
86
90template <typename TDomain, typename TAlgebra>
92(
93 const char* subsets
94)
95{
96 add_subset ("", subsets);
97}
98
102template <typename TDomain, typename TAlgebra>
104(
105 MathVector<dim> & value,
106 const char* function,
107 const char* subsets
108)
109{
110 add_subset (function, subsets);
111 m_vConstBCData.push_back (TConstBC (value, function, subsets));
112}
113
117template <typename TDomain, typename TAlgebra>
119(
120 std::vector<number> vValue,
121 const char* function,
122 const char* subsets
123)
124{
125 if (vValue.size () != (size_t) dim)
126 UG_THROW ("NedelecDirichletBC:"
127 "Wrong dimensionality of a Dirichlet BC. Specify " << dim << " values.");
128 MathVector<dim> value;
129 for (size_t i = 0; i < (size_t) dim; i++) value[i] = vValue[i];
130 add (value, function, subsets);
131}
132
136template <typename TDomain, typename TAlgebra>
138(
139 SmartPtr<UserData<MathVector<dim>, dim> > & func,
140 const char* function,
141 const char* subsets
142)
143{
144 add_subset (function, subsets);
145 m_vUserDataBCData.push_back (TUserDataBC (func, function, subsets));
146}
147
148#ifdef UG_FOR_LUA
152template <typename TDomain, typename TAlgebra>
154(
155 const char* name,
156 const char* function,
157 const char* subsets
158)
159{
160 if (! LuaUserData<MathVector<dim>, dim>::check_callback_returns (name))
161 UG_THROW ("NedelecDirichletBC: Illegal BC specification."
162 " A " << dim << "d vector-valued function should be specified");
163 SmartPtr<UserData<MathVector<dim>, dim> > sp_userData =
164 LuaUserDataFactory<MathVector<dim>, dim>::create (name);
165 add (sp_userData, function, subsets);
166}
167#endif
168
172template <typename TDomain, typename TAlgebra>
174(
175 SubsetGroup & dirichlet_ssgrp
176) const
177{
178 typedef std::map<std::string, std::vector<std::string> > t_ss_map;
179
180// Loop the subset names
181 t_ss_map::const_iterator iterEnd = m_mDirichletSS.end ();
182 for (t_ss_map::const_iterator iter = m_mDirichletSS.begin (); iter != iterEnd; ++iter)
183 dirichlet_ssgrp.add (iter->first);
184}
185
190template <typename TDomain, typename TAlgebra>
192(
193 FunctionGroup& functionGroup,
194 SubsetGroup& subsetGroup
195)
196{
197// only number of functions allowed
198 if (functionGroup.size () != 1)
199 UG_THROW("NedelecDirichletBC:"
200 " Only single functions function allowed in the specification of"
201 " Dirichlet values, but the following functions given:"
202 << functionGroup);
203
204// get subsethandler
205 ConstSmartPtr<ISubsetHandler> pSH = base_type::m_spApproxSpace->subset_handler ();
206
207// loop subsets
208 for(size_t si = 0; si < subsetGroup.size (); ++si)
209 {
210 // get the subset and function indices
211 const int subsetIndex = subsetGroup[si];
212 const size_t fct = functionGroup[0];
213
214 // check that subsetIndex is valid
215 if (subsetIndex < 0 || subsetIndex >= pSH->num_subsets ())
216 UG_THROW ("NedelecDirichletBC:"
217 " Invalid Subset Index " << subsetIndex << ". (Valid is"
218 " 0, .. , " << pSH->num_subsets() <<").");
219
220 // check if function exist
221 if (fct >= base_type::m_spApproxSpace->function_pattern()->num_fct ())
222 UG_THROW ("NedelecDirichletBC:"
223 " Function " << fct << " does not exist in pattern.");
224
225 // check that function is defined for segment
226 if (!base_type::m_spApproxSpace->function_pattern()->is_def_in_subset (fct, subsetIndex))
227 UG_THROW ("NedelecDirichletBC:"
228 " Function " << fct << " not defined on subset " << subsetIndex);
229 }
230}
231
236template <typename TDomain, typename TAlgebra>
237template <typename TUserData>
239(
240 std::map<int, std::vector<TUserData *> > & mvUserDataBndSegment,
241 std::vector<TUserData> & vUserData
242)
243{
244// Clear the extracted data
245 mvUserDataBndSegment.clear ();
246
247// Fill the map
248 for (size_t i = 0; i < vUserData.size (); ++i)
249 {
250 FunctionGroup fctGrp;
251
252 // create Function Group and Subset Group
253 try
254 {
255 vUserData[i].ssGrp = base_type::m_spApproxSpace->subset_grp_by_name (vUserData[i].ssName.c_str());
256 }
257 UG_CATCH_THROW(" Subsets '" << vUserData[i].ssName << "' not all contained in ApproximationSpace.");
258
259 try
260 {
261 fctGrp = base_type::m_spApproxSpace->fct_grp_by_name (vUserData[i].fctName.c_str());
262 }
263 UG_CATCH_THROW (" Functions '" << vUserData[i].fctName << "' not all contained in ApproximationSpace.");
264
265 // check functions and subsets
266 check_functions_and_subsets (fctGrp, vUserData[i].ssGrp);
267
268 // set the function
269 vUserData[i].fct = fctGrp[0];
270
271 // loop subsets
272 for (size_t si = 0; si < vUserData[i].ssGrp.size (); ++si)
273 mvUserDataBndSegment[vUserData[i].ssGrp[si]].push_back (&vUserData[i]);
274 }
275}
276
280template <typename TDomain, typename TAlgebra>
282{
283 typedef std::map<std::string, std::vector<std::string> > t_ss_map;
284
285// Loop the subset names
286 t_ss_map::const_iterator iterEnd = m_mDirichletSS.end ();
287 for (t_ss_map::const_iterator iter = m_mDirichletSS.begin (); iter != iterEnd; ++iter)
288 if (iter->second.size () != 0) // only if there are implicitly defined BCs
289 {
290 // Get the subset group
291 SubsetGroup ssGrp;
292 try
293 {
294 ssGrp = base_type::m_spApproxSpace->subset_grp_by_name (iter->first.c_str ());
295 }
296 UG_CATCH_THROW(" Subsets '" << iter->first << "' not all contained in ApproximationSpace.");
297
298 // Loop subsets
299 for (size_t i = 0; i < ssGrp.size (); i++)
300 {
301 FunctionGroup & f_grp = m_mZeroBC [ssGrp [i]];
302 f_grp.set_function_pattern (base_type::m_spApproxSpace->function_pattern ());
303 f_grp.add (iter->second);
304 }
305 }
306}
307
311template <typename TDomain, typename TAlgebra>
313(
314 matrix_type & J,
315 const vector_type & u,
317 int type,
318 number time,
320 const number s_a0
321)
322{
323 std::vector<DoFIndex> multInd;
324 SubsetGroup ssGrp (base_type::m_spApproxSpace->subset_handler ());
325 FunctionGroup fctGrp (base_type::m_spApproxSpace->function_pattern (), m_vDirichletFunc);
326
327// Get the Dirichlet subsets
329
330// Loop over the subsets
331 for (size_t i = 0; i < ssGrp.size (); i++)
332 {
333 int si = ssGrp [i];
334
335 // Loop the edges
336 t_edge_iterator iterEnd = dd->end<Edge> (si);
337 for (t_edge_iterator iter = dd->begin<Edge> (si); iter != iterEnd; iter++)
338 // Loop functions
339 for(size_t fct_i = 0; fct_i < fctGrp.size (); fct_i++)
340 {
341 // Get the multiindices
342 if (dd->inner_dof_indices (*iter, fctGrp [fct_i], multInd) != 1)
343 UG_THROW ("NedelecDirichletBC: "
344 "More than one DoF per edge. Not the Nedelec-type-1 element?");
345
346 // Set the Dirichlet row
347 SetDirichletRow (J, multInd[0]);
348 }
349 }
350}
351
355template <typename TDomain, typename TAlgebra>
357(
358 vector_type& d,
359 const vector_type& u,
361 int type,
362 number time,
364 const std::vector<number> * vScaleMass,
365 const std::vector<number> * vScaleStiff
366)
367{
368 std::vector<DoFIndex> multInd;
369 SubsetGroup ssGrp (base_type::m_spApproxSpace->subset_handler ());
370 FunctionGroup fctGrp (base_type::m_spApproxSpace->function_pattern (), m_vDirichletFunc);
371
372// Get the Dirichlet subsets
374
375// Loop over the subsets
376 for (size_t i = 0; i < ssGrp.size (); i++)
377 {
378 int si = ssGrp [i];
379
380 // Loop the edges
381 t_edge_iterator iterEnd = dd->end<Edge> (si);
382 for (t_edge_iterator iter = dd->begin<Edge> (si); iter != iterEnd; iter++)
383 // Loop functions
384 for(size_t fct_i = 0; fct_i < fctGrp.size (); fct_i++)
385 {
386 // Get the multiindices
387 if (dd->inner_dof_indices (*iter, fctGrp [fct_i], multInd) != 1)
388 UG_THROW ("NedelecDirichletBC: "
389 "More than one DoF per edge. Not the Nedelec-type-1 element?");
390
391 // Set the Dirichlet entry
392 DoFRef (d, multInd[0]) = 0.0;
393 }
394 }
395}
396
400template <typename TDomain, typename TAlgebra>
401template <typename TUserData>
403(
404 const std::vector<TUserData *> & vUserData,
405 int si,
406 vector_type& u,
408 number time
409)
410{
411 std::vector<DoFIndex> multInd;
412
413// Loop the edges
414 t_edge_iterator iterEnd = dd->end<Edge> (si);
415 for (t_edge_iterator iter = dd->begin<Edge> (si); iter != iterEnd; iter++)
416 // Loop dirichlet functions on this segment
417 for(size_t i = 0; i < vUserData.size (); ++i)
418 {
419 Edge * pEdge = *iter;
420 // Get the multiindices
421 if (dd->inner_dof_indices (pEdge, vUserData[i]->fct, multInd) != 1)
422 UG_THROW ("NedelecDirichletBC:"
423 "More than one DoF per edge. Not the Nedelec-type-1 element?");
424
425 // Set the Dirichlet entry
426 DoFRef (u, multInd[0]) = (* vUserData[i]) (pEdge, si, m_aaPos, time);
427 }
428}
429
434template <typename TDomain, typename TAlgebra>
436(
437 vector_type & u,
439 number time
440)
441{
442 typedef std::map<int, FunctionGroup> t_ss_map;
443
444 std::vector<DoFIndex> multInd;
445
446// Loop the subset names
447 t_ss_map::const_iterator ssIterEnd = m_mZeroBC.end ();
448 for (t_ss_map::const_iterator ssIter = m_mZeroBC.begin (); ssIter != ssIterEnd; ++ssIter)
449 {
450 const FunctionGroup & fctGrp = ssIter->second;
451 int si = ssIter->first;
452
453 // Loop the edges
454 t_edge_iterator iterEnd = dd->end<Edge> (si);
455 for (t_edge_iterator iter = dd->begin<Edge> (si); iter != iterEnd; iter++)
456 {
457 Edge * pEdge = *iter;
458
459 // Loop dirichlet functions on this segment
460 for(size_t i = 0; i < fctGrp.size (); ++i)
461 {
462 // Get the multiindices
463 if (dd->inner_dof_indices (pEdge, fctGrp[i], multInd) != 1)
464 UG_THROW ("NedelecDirichletBC:"
465 "More than one DoF per edge. Not the Nedelec-type-1 element?");
466
467 // Set the Dirichlet entry
468 DoFRef (u, multInd[0]) = 0;
469 }
470 }
471 }
472}
473
477template <typename TDomain, typename TAlgebra>
479(
480 vector_type & u,
482 int type,
483 number time
484)
485{
486 extract_data ();
487
488 // Loop boundary subsets (separately for every type of the BC specification)
489 {
490 typename std::map<int, std::vector<TConstBC*> >::const_iterator iter;
491 for (iter = m_mConstBC.begin (); iter != m_mConstBC.end (); ++iter)
492 adjust_solution<TConstBC> ((*iter).second, (*iter).first, u, dd, time);
493 }
494 {
495 typename std::map<int, std::vector<TUserDataBC*> >::const_iterator iter;
496 for (iter = m_mUserDataBC.begin (); iter != m_mUserDataBC.end (); ++iter)
497 adjust_solution<TUserDataBC> ((*iter).second, (*iter).first, u, dd, time);
498 }
499
500 adjust_solution_implicit (u, dd, time);
501}
502
503} // end namespace Electromagnetism
504} // end namespace ug
505
506/* End of File */
void adjust_jacobian(matrix_type &J, const vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol=SPNULL, const number s_a0=1.0)
sets a unity row for all dirichlet indices
Definition nedelec_dirichlet_impl.h:313
algebra_type::matrix_type matrix_type
type of algebra matrix
Definition nedelec_dirichlet.h:107
algebra_type::vector_type vector_type
type of algebra vector
Definition nedelec_dirichlet.h:110
void add_subset(const char *f_name, const char *ss_names)
inserts a subset to the list of the Dirichlet subsets
Definition nedelec_dirichlet_impl.h:50
void extract_data()
Definition nedelec_dirichlet.h:294
void adjust_solution(const std::vector< TUserData * > &vUserData, int si, vector_type &u, ConstSmartPtr< DoFDistribution > dd, number time)
Sets the specified Dirichlet values in the solution.
Definition nedelec_dirichlet_impl.h:403
void add(MathVector< dim > &value, const char *function, const char *subsets)
adds a constant Dirichlet value on a subset
Definition nedelec_dirichlet_impl.h:104
void adjust_solution_implicit(vector_type &u, ConstSmartPtr< DoFDistribution > dd, number time)
Sets the zero Dirichlet values for the implicitly specified BC.
Definition nedelec_dirichlet_impl.h:436
void extract_implicit()
Composes the list of the implicit BC.
Definition nedelec_dirichlet_impl.h:281
void check_functions_and_subsets(FunctionGroup &functionGroup, SubsetGroup &subsetGroup)
Verifies whether the string input represents legal data.
Definition nedelec_dirichlet_impl.h:192
void adjust_defect(vector_type &d, const vector_type &u, ConstSmartPtr< DoFDistribution > dd, int type, number time=0.0, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol=SPNULL, const std::vector< number > *vScaleMass=NULL, const std::vector< number > *vScaleStiff=NULL)
sets a zero value in the defect for all dirichlet indices
Definition nedelec_dirichlet_impl.h:357
void add_0(const char *subsets)
adds a zero Dirichlet value on a subset
Definition nedelec_dirichlet_impl.h:92
virtual void get_dirichlet_subsets(SubsetGroup &dirichlet_ssgrp) const
composes an array of subset indices of the Dirichlet boundary
Definition nedelec_dirichlet_impl.h:174
size_t size() const
void add(const char *name)
void set_function_pattern(ConstSmartPtr< FunctionPattern > spFuncPattern)
void add(const char *name)
size_t size() const
SmartPtr< TSubsetHandler > subset_handler()
Variant::Type type()
vector< string > TokenizeString(const char *str, const char delimiter=',')
#define UG_CATCH_THROW(msg)
#define UG_THROW(msg)
double number
const number & DoFRef(const TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
void SetDirichletRow(TMatrix &mat, const DoFIndex &ind)
structure for the Dirichlet BC that is constant over a patch
Definition nedelec_dirichlet.h:165
Structure for the Dirichlet BC that are given by a function.
Definition nedelec_dirichlet.h:207