Plugins
eddy_current_cmd_impl.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 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 and the commands for the eddy current simulations.
35  */
36 
37 // basic ug4 headers
38 #include "common/common.h"
39 #ifdef UG_FOR_LUA
41 #endif
42 
43 /* Discretization's headers: */
44 #include "../nedelec_local_ass.h"
45 
46 namespace ug{
47 namespace Electromagnetism{
48 
49 //---- Computation of the power of the electromagnetic field ----//
50 
52 
62 template <typename TGridFunc, typename TElem>
64 {
65 private:
66  typedef typename TGridFunc::domain_type domain_type;
69 
73 
75  inline static void calc_elem_power
76  (
77  TGridFunc * pJGGF,
78  size_t JG_fct[],
79  TGridFunc * pEGF,
80  size_t E_fct[],
81  position_accessor_type & aaPos,
82  TElem * pElem,
83  number pow[]
84  )
85  {
86  // Get coordinates of the corners
87  position_type corners [numCorners];
88  for (size_t co = 0; co < numCorners; co++)
89  corners [co] = aaPos [pElem->vertex (co)];
90 
91  // Assemble the local mass matrix
92  number M [maxEdges] [maxEdges];
94  (pEGF->domain().get (), pElem, corners, M);
95 
96  // Get the components of the grid functions
97  number Re_JG_val [numEdges], Im_JG_val [numEdges];
98  number Re_E_val [numEdges], Im_E_val [numEdges];
99  std::vector<DoFIndex> ind;
100 
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]);
105 
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]);
110 
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]);
115 
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]);
120 
121  // Compute the integral
122  number Re_ME [numEdges], Im_ME [numEdges];
123 
124  memset (Re_ME, 0, numEdges * sizeof (number));
125  memset (Im_ME, 0, numEdges * sizeof (number));
126  for (size_t i = 0; i < numEdges; i++)
127  for (size_t j = 0; j < numEdges; j++)
128  {
129  Re_ME [i] += M [i] [j] * Re_E_val [j];
130  Im_ME [i] += M [i] [j] * Im_E_val [j];
131  }
132 
133  // a) the real part
134  number Re_pow = 0;
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];
137  // b) the imaginary part
138  number Im_pow = 0;
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];
141 
142  // Done
143  pow[0] += Re_pow; pow[1] += Im_pow;
144 
145  //-- For debugging only
146  // UG_LOG ("==> subset " << pEGF->domain()->subset_handler()->get_subset_index (pElem)
147  // << ": Re = " << Re_pow << ", Im = " << Im_pow << "\n");
148  //--
149  }
150 
151 public:
152 
154  static void calc_power
155  (
156  TGridFunc * pJGGF,
157  size_t JG_fct[],
158  SubsetGroup & JG_ssg,
159  TGridFunc * pEGF,
160  size_t E_fct[],
161  number pow[]
162  )
163  {
164  typedef typename TGridFunc::template traits<TElem>::const_iterator t_elem_iterator;
165 
166  position_accessor_type & aaPos = pEGF->domain()->position_accessor();
167 
168  // loop the subdomains
169  for (size_t i = 0; i < JG_ssg.size (); i++)
170  {
171  int si = JG_ssg [i];
172 
173  // loop all the elements of the given type in the subset
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);
179  }
180  }
181 };
182 
184 template <typename TGridFunc>
186 {
187  TGridFunc * m_pJGGF;
188  size_t m_JG_fct[2];
190  TGridFunc * m_pEGF;
191  size_t m_E_fct[2];
193 
194 public:
195 
198  (
199  TGridFunc * pJGGF,
200  size_t JG_fct[],
201  SubsetGroup & JG_ssg,
202  TGridFunc * pEGF,
203  size_t E_fct[],
204  number pow[]
205  )
206  : m_pJGGF (pJGGF), m_JG_ssg (JG_ssg), m_pEGF (pEGF), m_pow (pow)
207  {
208  m_JG_fct[0] = JG_fct[0]; m_JG_fct[1] = JG_fct[1];
209  m_E_fct[0] = E_fct[0]; m_E_fct[1] = E_fct[1];
210  m_pow[0] = m_pow[1] = 0;
211  }
212 
213  template <typename TElem> void operator() (TElem &)
214  {
217  }
218 };
219 
221 template <typename TGridFunc>
223 (
224  TGridFunc * pJGGF,
225  size_t JG_fct[],
226  SubsetGroup & JG_ssg,
227  TGridFunc * pEGF,
228  size_t E_fct[],
229  number pow[]
230 )
231 {
232  typedef typename TGridFunc::domain_type domain_type;
233  static const int dim = domain_type::dim;
234  typedef typename domain_traits<dim>::DimElemList ElemList;
235 
236  boost::mpl::for_each<ElemList>
237  (CalcVolPowerHelperClass<TGridFunc> (pJGGF, JG_fct, JG_ssg, pEGF, E_fct, pow));
238 
239  pow [0] /= -2;
240  pow [1] /= -2;
241 }
242 
244 template <typename TGridFunc>
246 (
247  SmartPtr<TGridFunc> spJGGF,
248  const char* JG_cmps,
249  const char* JG_ss,
250  SmartPtr<TGridFunc> spEGF,
251  const char* E_cmps
252 )
253 {
254 // Get functions by names
255  size_t JG_fcts [2];
256  {
257  std::string fctString = std::string (JG_cmps);
258  std::vector<std::string> tokens;
259  TokenizeString (fctString, 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++)
265  RemoveWhitespaceFromString (tokens [i]);
266 
267  for (int i = 0; i < 2; i++)
268  try
269  {
270  JG_fcts [i] = spJGGF->fct_id_by_name (tokens[i].c_str ());
271  }
272  UG_CATCH_THROW ("EddyCurrent: Cannot find symbolic function "
273  "component for the name '" << tokens[i] << "' (in JG).");
274  }
275  size_t E_fcts [2];
276  {
277  std::string fctString = std::string (E_cmps);
278  std::vector<std::string> tokens;
279  TokenizeString (fctString, 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++)
285  RemoveWhitespaceFromString (tokens [i]);
286 
287  for (int i = 0; i < 2; i++)
288  try
289  {
290  E_fcts [i] = spEGF->fct_id_by_name (tokens[i].c_str ());
291  }
292  UG_CATCH_THROW ("EddyCurrent: Cannot find symbolic function "
293  "component for the name '" << tokens[i] << "' (in E).");
294  }
295 
296  SubsetGroup JG_ssg;
297  {
298  std::string ssString = std::string (JG_ss);
299  std::vector<std::string> tokens;
300  TokenizeString (ssString, tokens, ',');
301  for (size_t i = 0; i < tokens.size (); i++)
302  RemoveWhitespaceFromString (tokens [i]);
303  JG_ssg.set_subset_handler (spJGGF->subset_handler ());
304  JG_ssg.add (tokens);
305  }
306 
307 // Compute the power
308  number power [2];
309  calc_power (spJGGF.get (), JG_fcts, JG_ssg, spEGF.get (), E_fcts, power);
310 
311 // Print the result
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");
315 }
316 
317 //---- Computation of magnetic flux through a cross-section ----//
318 
320 
327 template <typename TGridFunc>
329 (
330  const TGridFunc * gfE,
331  const size_t fct [2],
332  const SubsetGroup & ss_grp,
333  const MathVector<TGridFunc::dim> & Normal,
334  const typename TGridFunc::domain_type::position_type & base_pnt,
335  const size_t n_pnt,
336  const typename TGridFunc::domain_type::position_type & d_pnt,
337  number emf [2]
338 )
339 {
340  typedef typename TGridFunc::domain_type domain_type;
341  typedef typename domain_type::position_type position_type;
342  static const int dim = TGridFunc::dim;
343  typedef typename domain_traits<dim>::grid_base_object elem_type;
344  typedef typename TGridFunc::template traits<elem_type>::const_iterator t_elem_iterator;
345  static const size_t maxCorners = domain_traits<dim>::MaxNumVerticesOfElem;
346  static const size_t maxEdges = (size_t) element_list_traits<typename domain_traits<dim>::DimElemList>::maxEdges;
347 
348  position_type corners [maxCorners];
349  std::vector<DoFIndex> ind;
350  number dofValues [2] [maxEdges];
351 
352 // check the function space of the grid function
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)");
356 
357 // check whether the functions are defined in the specified subsets
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.");
362 
363 // normalize the normal
364  number normal_len = VecLength (Normal);
365  if (normal_len < 1e-16)
366  UG_THROW ("EddyCurrent: Illegal normal vector specified (almost zero).");
367  MathVector<dim> normal (Normal);
368  normal *= 1 / normal_len;
369 
370 // get the domain
371  const typename TGridFunc::domain_type * domain = gfE->domain().get ();
372 
373 // get position accessor
374  const typename domain_type::position_accessor_type& aaPos = domain->position_accessor ();
375 
376 // compute the flux for every subset
377  number total_area = 0;
378  emf[0] = emf[1] = 0;
379  for (size_t ssgi = 0; ssgi < ss_grp.size (); ssgi++)
380  {
381  int si = ss_grp [ssgi];
382 
383  // Loop the (full-dim) elements
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++)
386  {
387  elem_type * elem = *iter;
388 
389  // get the corner coordinates
390  for (size_t co = 0; co < elem->num_vertices (); co++)
391  corners [co] = aaPos [elem->vertex (co)];
392 
393  // get the DoFs
394  for (size_t part = 0; part < 2; part++)
395  {
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]);
400  }
401 
402  // loop over the windings
403  position_type pnt = base_pnt;
404  for (size_t i = 0; i < n_pnt; i++, pnt += d_pnt)
405  {
406  number elem_flux [2];
407  number elem_area;
409  (domain, (GridObject *) elem, corners, dofValues[0], normal, pnt, elem_flux[0])) != 0.0)
410  {
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;
416  }
417  }
418  }
419  }
420 
421 // done
422  return total_area;
423 };
424 
426 template <typename TGridFunc>
428 (
429  SmartPtr<TGridFunc> spGF,
430  const char* cmps,
431  const char* subsets,
432  const std::vector<number>& Normal,
433  const std::vector<number>& base_pnt,
434  const size_t n_pnt,
435  const std::vector<number>& d_pnt
436 )
437 {
438  typedef typename TGridFunc::domain_type domain_type;
439  typedef typename domain_type::position_type position_type;
440  static const int dim = TGridFunc::dim;
441 
442 // check the sizes of the arrays and copy the values
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.");
445  MathVector<dim> normal_vec;
446  position_type base_pnt_coord, d_pnt_coord;
447  for (int i = 0; i < dim; i++)
448  {
449  normal_vec[i] = Normal[i];
450  base_pnt_coord[i] = base_pnt[i];
451  d_pnt_coord[i] = d_pnt[i];
452  }
453 
454 // create the subset group
455  if (subsets == NULL)
456  UG_THROW ("EddyCurrent: No subsets specified.");
457  SubsetGroup ssGrp (spGF->domain()->subset_handler ());
458  {
459  std::string ssString = std::string (subsets);
460  std::vector<std::string> tokens;
461  TokenizeString (ssString, tokens, ',');
462  for (size_t i = 0; i < tokens.size (); i++)
463  RemoveWhitespaceFromString (tokens [i]);
464  ssGrp.add (tokens);
465  }
466 
467 // get functions by names
468  size_t fcts [2];
469  {
470  std::string fctString = std::string (cmps);
471  std::vector<std::string> tokens;
472  TokenizeString (fctString, 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++)
478  RemoveWhitespaceFromString (tokens [i]);
479 
480  for (int i = 0; i < 2; i++)
481  try
482  {
483  fcts [i] = spGF->fct_id_by_name (tokens[i].c_str ());
484  }
485  UG_CATCH_THROW ("EddyCurrent: Cannot find symbolic function "
486  "component for the name '" << tokens[i] << "'.");
487  }
488 
489 // compute the flux and the area
490  number emf [2], area;
491  area = calc_EMF (spGF.get(), fcts, ssGrp, normal_vec,
492  base_pnt_coord, n_pnt, d_pnt_coord, emf);
493 
494 // print the result
495  UG_LOG ("--> Electromotive force in the coil with " << n_pnt << " windings (total area "
496  << area << "): " << emf[0] << " + " << emf[1] << " I.\n");
497 };
498 
499 } // end namespace Electromagnetism
500 } // end namespace ug
501 
502 /* End of File */
T * get()
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)
size_t size() const
void set_subset_handler(ConstSmartPtr< ISubsetHandler > sh)
Grid::VertexAttachmentAccessor< position_attachment_type > position_accessor_type
MathVector< dim > position_type
static const int dim
normal
vector< string > TokenizeString(const char *str, const char delimiter=',')
void RemoveWhitespaceFromString(std::string &string)
#define UG_CATCH_THROW(msg)
#define UG_THROW(msg)
#define UG_LOG(msg)
double number
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