44 #include "../nedelec_local_ass.h"
47 namespace Electromagnetism{
62 template <
typename TGr
idFunc,
typename TElem>
89 corners [co] = aaPos [pElem->vertex (co)];
94 (pEGF->domain().get (), pElem, corners, M);
99 std::vector<DoFIndex> ind;
101 if (pJGGF->dof_indices (pElem, JG_fct[0], ind) !=
numEdges)
102 UG_THROW (
"J_G grid function type mismatch.");
103 for (
size_t dof = 0; dof <
numEdges; dof++)
104 Re_JG_val[dof] =
DoFRef (*pJGGF, ind[dof]);
106 if (pJGGF->dof_indices (pElem, JG_fct[1], ind) !=
numEdges)
107 UG_THROW (
"J_G grid function type mismatch.");
108 for (
size_t dof = 0; dof <
numEdges; dof++)
109 Im_JG_val[dof] =
DoFRef (*pJGGF, ind[dof]);
111 if (pEGF->dof_indices (pElem, E_fct[0], ind) !=
numEdges)
112 UG_THROW (
"E grid function type mismatch.");
113 for (
size_t dof = 0; dof <
numEdges; dof++)
114 Re_E_val[dof] =
DoFRef (*pEGF, ind[dof]);
116 if (pEGF->dof_indices (pElem, E_fct[1], ind) !=
numEdges)
117 UG_THROW (
"E grid function type mismatch.");
118 for (
size_t dof = 0; dof <
numEdges; dof++)
119 Im_E_val[dof] =
DoFRef (*pEGF, ind[dof]);
126 for (
size_t i = 0; i <
numEdges; i++)
127 for (
size_t j = 0; j <
numEdges; j++)
129 Re_ME [i] += M [i] [j] * Re_E_val [j];
130 Im_ME [i] += M [i] [j] * Im_E_val [j];
135 for (
size_t i = 0; i <
numEdges; i++)
136 Re_pow += Re_JG_val[i] * Re_ME[i] + Im_JG_val[i] * Im_ME[i];
139 for (
size_t i = 0; i <
numEdges; i++)
140 Im_pow += Re_JG_val[i] * Im_ME[i] - Re_ME[i] * Im_JG_val[i];
143 pow[0] += Re_pow; pow[1] += Im_pow;
169 for (
size_t i = 0; i < JG_ssg.
size (); i++)
174 t_elem_iterator elem_iter = pEGF->template begin<TElem> (si);
175 t_elem_iterator end_iter = pEGF->template end<TElem> (si);
176 for (; elem_iter != end_iter; ++elem_iter)
178 (pJGGF, JG_fct, pEGF, E_fct, aaPos, *elem_iter, pow);
184 template <
typename TGr
idFunc>
221 template <
typename TGr
idFunc>
232 typedef typename TGridFunc::domain_type domain_type;
236 boost::mpl::for_each<ElemList>
244 template <
typename TGr
idFunc>
257 std::string fctString = std::string (JG_cmps);
258 std::vector<std::string> tokens;
260 if (tokens.size () != 2)
261 UG_THROW(
"EddyCurrent: Needed 2 components "
262 "in symbolic function names (for Re and Im) for JG, "
263 "but given: " << JG_cmps);
264 for (
size_t i = 0; i < tokens.size (); i++)
267 for (
int i = 0; i < 2; i++)
270 JG_fcts [i] = spJGGF->fct_id_by_name (tokens[i].c_str ());
273 "component for the name '" << tokens[i] <<
"' (in JG).");
277 std::string fctString = std::string (E_cmps);
278 std::vector<std::string> tokens;
280 if (tokens.size () != 2)
281 UG_THROW(
"EddyCurrent: Needed 2 components "
282 "in symbolic function names (for Re and Im) for E, "
283 "but given: " << E_cmps);
284 for (
size_t i = 0; i < tokens.size (); i++)
287 for (
int i = 0; i < 2; i++)
290 E_fcts [i] = spEGF->fct_id_by_name (tokens[i].c_str ());
293 "component for the name '" << tokens[i] <<
"' (in E).");
298 std::string ssString = std::string (JG_ss);
299 std::vector<std::string> tokens;
301 for (
size_t i = 0; i < tokens.size (); i++)
312 UG_LOG (
"--> Power of the electromagnetic field: "
313 << power [0] <<
" + " << power [1] <<
" I.\n");
314 UG_LOG (
"Note: Contribution of the boundary is not taken into account.\n");
327 template <
typename TGr
idFunc>
330 const TGridFunc * gfE,
331 const size_t fct [2],
340 typedef typename TGridFunc::domain_type domain_type;
344 typedef typename TGridFunc::template traits<elem_type>::const_iterator t_elem_iterator;
349 std::vector<DoFIndex> ind;
350 number dofValues [2] [maxEdges];
353 for (
size_t part = 0; part < 2; part++)
354 if (gfE->local_finite_element_id(fct[part]).type () !=
LFEID::NEDELEC)
355 UG_THROW (
"EddyCurrent: Illegal discretization space of the grid func. (must be NEDELEC)");
358 for (
size_t ssgi = 0; ssgi < ss_grp.
size (); ssgi++)
359 for (
size_t i = 0; i < 2; i++)
360 if (! gfE->is_def_in_subset (i, ss_grp [ssgi]))
361 UG_THROW (
"EddyCurrent: Grid function is undefined in some of the subsets.");
365 if (normal_len < 1e-16)
366 UG_THROW (
"EddyCurrent: Illegal normal vector specified (almost zero).");
371 const typename TGridFunc::domain_type * domain = gfE->domain().get ();
379 for (
size_t ssgi = 0; ssgi < ss_grp.
size (); ssgi++)
381 int si = ss_grp [ssgi];
384 t_elem_iterator iterEnd = gfE->template end<elem_type> (si);
385 for (t_elem_iterator iter = gfE->template begin<elem_type> (si); iter != iterEnd; iter++)
387 elem_type * elem = *iter;
390 for (
size_t co = 0; co < elem->num_vertices (); co++)
391 corners [co] = aaPos [elem->vertex (co)];
394 for (
size_t part = 0; part < 2; part++)
396 if (gfE->dof_indices (elem, fct[part], ind) > maxEdges)
397 UG_THROW (
"EddyCurrent: Illegal number of DoFs per element.");
398 for (
size_t sh = 0; sh < ind.size (); sh++)
399 dofValues [part] [sh] =
DoFRef (*gfE, ind [sh]);
404 for (
size_t i = 0; i < n_pnt; i++, pnt += d_pnt)
409 (domain, (
GridObject *) elem, corners, dofValues[0],
normal, pnt, elem_flux[0])) != 0.0)
412 (domain, (
GridObject *) elem, corners, dofValues[1],
normal, pnt, elem_flux[1]);
413 emf[0] += elem_flux[0];
414 emf[1] += elem_flux[1];
415 total_area += elem_area;
426 template <
typename TGr
idFunc>
432 const std::vector<number>& Normal,
433 const std::vector<number>& base_pnt,
435 const std::vector<number>& d_pnt
438 typedef typename TGridFunc::domain_type domain_type;
443 if (Normal.size () != (
size_t)
dim || base_pnt.size () != (
size_t)
dim || d_pnt.size () != (
size_t)
dim)
444 UG_THROW (
"EddyCurrent: Illegal sizes of the vectors in the arguments.");
447 for (
int i = 0; i <
dim; i++)
449 normal_vec[i] = Normal[i];
450 base_pnt_coord[i] = base_pnt[i];
451 d_pnt_coord[i] = d_pnt[i];
456 UG_THROW (
"EddyCurrent: No subsets specified.");
457 SubsetGroup ssGrp (spGF->domain()->subset_handler ());
459 std::string ssString = std::string (subsets);
460 std::vector<std::string> tokens;
462 for (
size_t i = 0; i < tokens.size (); i++)
470 std::string fctString = std::string (cmps);
471 std::vector<std::string> tokens;
473 if (tokens.size () != 2)
474 UG_THROW(
"EddyCurrent: Needed 2 components "
475 "in symbolic function names (for Re and Im), "
476 "but given: " << cmps);
477 for (
size_t i = 0; i < tokens.size (); i++)
480 for (
int i = 0; i < 2; i++)
483 fcts [i] = spGF->fct_id_by_name (tokens[i].c_str ());
486 "component for the name '" << tokens[i] <<
"'.");
491 area =
calc_EMF (spGF.
get(), fcts, ssGrp, normal_vec,
492 base_pnt_coord, n_pnt, d_pnt_coord, emf);
495 UG_LOG (
"--> Electromotive force in the coil with " << n_pnt <<
" windings (total area "
496 << area <<
"): " << emf[0] <<
" + " << emf[1] <<
" I.\n");
Helper class for the integration of the power.
Definition: eddy_current_cmd_impl.h:64
domain_type::position_type position_type
Definition: eddy_current_cmd_impl.h:68
TGridFunc::domain_type domain_type
Definition: eddy_current_cmd_impl.h:66
domain_type::position_accessor_type position_accessor_type
Definition: eddy_current_cmd_impl.h:67
static const size_t numEdges
Definition: eddy_current_cmd_impl.h:71
static const size_t maxEdges
Definition: eddy_current_cmd_impl.h:72
static const size_t numCorners
Definition: eddy_current_cmd_impl.h:70
static void calc_power(TGridFunc *pJGGF, size_t JG_fct[], SubsetGroup &JG_ssg, TGridFunc *pEGF, size_t E_fct[], number pow[])
Computes the integral for all elements of one type.
Definition: eddy_current_cmd_impl.h:155
static void calc_elem_power(TGridFunc *pJGGF, size_t JG_fct[], TGridFunc *pEGF, size_t E_fct[], position_accessor_type &aaPos, TElem *pElem, number pow[])
computes the integral over one element
Definition: eddy_current_cmd_impl.h:76
Helper class for the computation of the power of the electromagnetic field.
Definition: eddy_current_cmd_impl.h:186
CalcVolPowerHelperClass(TGridFunc *pJGGF, size_t JG_fct[], SubsetGroup &JG_ssg, TGridFunc *pEGF, size_t E_fct[], number pow[])
class constructor
Definition: eddy_current_cmd_impl.h:198
TGridFunc * m_pJGGF
Definition: eddy_current_cmd_impl.h:187
number * m_pow
Definition: eddy_current_cmd_impl.h:192
SubsetGroup & m_JG_ssg
Definition: eddy_current_cmd_impl.h:189
void operator()(TElem &)
Definition: eddy_current_cmd_impl.h:213
size_t m_JG_fct[2]
Definition: eddy_current_cmd_impl.h:188
TGridFunc * m_pEGF
Definition: eddy_current_cmd_impl.h:190
size_t m_E_fct[2]
Definition: eddy_current_cmd_impl.h:191
Interpolation of the Nedelec dofs and their curls.
Definition: nedelec_local_ass.h:399
static number curl_flux(const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< WDim > &normal, const position_type pnt, number &flux)
computes flux of the curl through a given plane in an element
Definition: nedelec_local_ass.h:457
Tool kit for the Whitney-1 (Nedelec) based FE discretization of the rot-rot operators.
Definition: nedelec_local_ass.h:85
static void local_mass(const TDomain *domain, TElem *elem, const position_type *corners, number M[maxNumEdges][maxNumEdges])
assembles the local mass matrix of the Nedelec elements
Definition: nedelec_local_ass.h:132
void add(const char *name)
void set_subset_handler(ConstSmartPtr< ISubsetHandler > sh)
Grid::VertexAttachmentAccessor< position_attachment_type > position_accessor_type
MathVector< dim > position_type
vector< string > TokenizeString(const char *str, const char delimiter=',')
void RemoveWhitespaceFromString(std::string &string)
#define UG_CATCH_THROW(msg)
vector_t::value_type VecLength(const vector_t &v)
void calc_power(TGridFunc *pJGGF, size_t JG_fct[], SubsetGroup &JG_ssg, TGridFunc *pEGF, size_t E_fct[], number pow[])
Computes the power of the electromagnetic field (up to the contribution of the boundary)
Definition: eddy_current_cmd_impl.h:223
void CalcEMF(SmartPtr< TGridFunc > spGF, const char *cmps, const char *subsets, const std::vector< number > &Normal, const std::vector< number > &base_pnt, const size_t n_pnt, const std::vector< number > &d_pnt)
Prints of the electromotive force due to the time variation of the magnetic flux through windings of ...
Definition: eddy_current_cmd_impl.h:428
void CalcPower(SmartPtr< TGridFunc > spJGGF, const char *JG_cmps, const char *JG_ss, SmartPtr< TGridFunc > spEGF, const char *E_cmps)
Prints the (complex-valued) power of the electromagnetic field.
Definition: eddy_current_cmd_impl.h:246
number calc_EMF(const TGridFunc *gfE, const size_t fct[2], const SubsetGroup &ss_grp, const MathVector< TGridFunc::dim > &Normal, const typename TGridFunc::domain_type::position_type &base_pnt, const size_t n_pnt, const typename TGridFunc::domain_type::position_type &d_pnt, number emf[2])
Computation of the magnetic flux through windings of a coil.
Definition: eddy_current_cmd_impl.h:329
const number & DoFRef(const TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
geometry_traits< TElem >::const_iterator const_iterator