Plugins
ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra > Class Template Reference

FE-discretization of the time-harmonic E-based formulation of the eddy current model. More...

#include <eddy_current_e_nedelec.h>

+ Inheritance diagram for ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >:

Classes

struct  RegisterLocalDiscr
 
struct  tGeneratorCurrent
 class for a generator current (source) in a subdomain More...
 

Public Member Functions

 EddyCurrent_E_Nedelec (const char *functions, ConstSmartPtr< EMaterial< domain_type > > spSubsetData, number frequency)
 class constructor More...
 
void set_generator_current (SmartPtr< TGridFunction > spgfJG, const char *cmp)
 adds a generator current item \( \mathbf{J}_{G,h} \) to the discretization More...
 
void set_generator_current (SmartPtr< TGridFunction > spgfJG, const char *cmp, const char *ss_names=NULL)
 adds a generator current item \( \mathbf{J}_{G,h} \) to the discretization More...
 
- Public Member Functions inherited from ug::IElemDisc< TDomain >
void add_elem_modifier (SmartPtr< IElemDiscModifier< TDomain > > elemModifier)
 
std::vector< SmartPtr< IElemDiscModifier< TDomain > > > & get_elem_modifier ()
 
 IElemDisc (const char *functions, const char *subsets)
 
 IElemDisc (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset)
 
- Public Member Functions inherited from IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >
virtual void add_def_A_elem (LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
virtual void add_def_A_expl_elem (LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
virtual void add_def_M_elem (LocalVector &d, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
virtual void add_jac_A_elem (LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
virtual void add_jac_M_elem (LocalMatrix &J, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
virtual void add_rhs_elem (LocalVector &rhs, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
IElemDisc< TDomain > & asLeaf ()
 
void check_roid (ReferenceObjectID roid, int discType)
 
void do_add_def_A_elem (LocalVector &d, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void do_add_def_A_expl_elem (LocalVector &d, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void do_add_def_M_elem (LocalVector &d, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void do_add_jac_A_elem (LocalMatrix &J, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void do_add_jac_M_elem (LocalMatrix &J, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void do_add_rhs_elem (LocalVector &rhs, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void do_fsh_elem_loop ()
 
void do_fsh_timestep (const number time, VectorProxyBase *u, size_t algebra_id)
 
void do_fsh_timestep_elem (const number time, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void do_prep_elem (LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
 
void do_prep_elem_loop (const ReferenceObjectID roid, const int si)
 
void do_prep_timestep (number future_time, const number time, VectorProxyBase *u, size_t algebra_id)
 
void do_prep_timestep_elem (const number time, LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
virtual void fsh_elem_loop ()
 
virtual void fsh_timestep (number time, VectorProxyBase *u)
 
virtual void fsh_timestep_elem (const number time, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
 IElemAssembleFuncs ()
 
virtual void post_assemble_loop ()
 
virtual void prep_assemble_loop ()
 
virtual void prep_elem (const LocalVector &u, GridObject *elem, const ReferenceObjectID roid, const MathVector< dim > vCornerCoords[])
 
virtual void prep_elem_loop (const ReferenceObjectID roid, const int si)
 
virtual void prep_timestep (number future_time, number time, VectorProxyBase *u)
 
virtual void prep_timestep_elem (const number time, const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void set_roid (ReferenceObjectID id, int discType)
 
virtual ~IElemAssembleFuncs ()
 
- Public Member Functions inherited from ug::IElemError< class >
 IElemError (const char *functions, const char *subsets)
 
 IElemError (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset)
 
- Public Member Functions inherited from ug::IElemDiscBase< class >
SmartPtr< ApproximationSpace< TDomain > > approx_space ()
 
ConstSmartPtr< ApproximationSpace< TDomain > > approx_space () const
 
void check_setup (bool bNonRegularGrid)
 
void clear_imports ()
 
SmartPtr< TDomain > domain ()
 
ConstSmartPtr< TDomain > domain () const
 
const FunctionGroupfunction_group () const
 
ConstSmartPtr< FunctionPatternfunction_pattern () const
 
IDataImport< dim > & get_import (size_t i)
 
 IElemDiscBase (const char *functions="", const char *subsets="")
 
 IElemDiscBase (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset)
 
bool is_time_dependent () const
 
bool local_time_series_needed ()
 
const LocalVectorTimeSerieslocal_time_solutions () const
 
const FunctionIndexMappingmap () const
 
number mass_scale () const
 
number mass_scale () const
 
number mass_scale (const size_t timePoint) const
 
number mass_scale (const size_t timePoint) const
 
const std::vector< number > & mass_scales () const
 
const std::vector< number > & mass_scales () const
 
size_t num_fct () const
 
size_t num_imports () const
 
size_t num_subsets () const
 
void register_import (IDataImport< dim > &Imp)
 
virtual bool requests_local_time_series ()
 
void set_approximation_space (SmartPtr< ApproximationSpace< TDomain > > approxSpace)
 
void set_functions (const std::string &functions)
 
void set_functions (const std::vector< std::string > &functions)
 
void set_stationary ()
 
void set_stationary (bool bStationaryForced=true)
 
void set_subsets (const std::string &subsets)
 
void set_subsets (const std::vector< std::string > &subsets)
 
void set_time_dependent (LocalVectorTimeSeries &locTimeSeries, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff)
 
void set_time_independent ()
 
void set_time_point (const size_t timePoint)
 
number stiff_scale () const
 
number stiff_scale () const
 
number stiff_scale (const size_t timePoint) const
 
number stiff_scale (const size_t timePoint) const
 
const std::vector< number > & stiff_scales () const
 
const std::vector< number > & stiff_scales () const
 
TDomain::subset_handler_typesubset_handler ()
 
const TDomain::subset_handler_typesubset_handler () const
 
const std::vector< std::string > & symb_fcts () const
 
const std::vector< std::string > & symb_subsets () const
 
number time () const
 
size_t time_point () const
 
virtual int type () const
 
virtual bool use_hanging () const
 
virtual ~IElemDiscBase ()
 
- Public Member Functions inherited from IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain >
IElemDisc< TDomain > & asLeaf ()
 
void check_roid (ReferenceObjectID roid, int discType)
 
virtual void compute_err_est_A_elem (const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
 
virtual void compute_err_est_M_elem (const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
 
virtual void compute_err_est_rhs_elem (GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
 
void do_compute_err_est_A_elem (LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
 
void do_compute_err_est_M_elem (LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
 
void do_compute_err_est_rhs_elem (GridObject *elem, const MathVector< dim > vCornerCoords[], const number &scale)
 
void do_fsh_err_est_elem_loop ()
 
void do_prep_err_est_elem (LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
void do_prep_err_est_elem_loop (const ReferenceObjectID roid, const int si)
 
virtual SmartPtr< IErrEstData< TDomain > > err_est_data ()
 
bool err_est_enabled () const
 
virtual void fsh_err_est_elem_loop ()
 
 IElemEstimatorFuncs ()
 
virtual void prep_err_est_elem (const LocalVector &u, GridObject *elem, const MathVector< dim > vCornerCoords[])
 
virtual void prep_err_est_elem_loop (const ReferenceObjectID roid, const int si)
 
void set_error_estimator (SmartPtr< IErrEstData< TDomain > > ee)
 
void set_roid (ReferenceObjectID id, int discType)
 
virtual ~IElemEstimatorFuncs ()
 

Private Types

typedef IElemDisc< TDomain > base_type
 base class type More...
 
typedef base_type::domain_type domain_type
 domain type More...
 
typedef base_type::position_type position_type
 position type More...
 
typedef GridFunction< TDomain, TAlgebra > TGridFunction
 type of grid functions (used for the sources) More...
 
typedef EddyCurrent_E_Nedelec< TDomain, TAlgebra > this_type
 own type More...
 

Private Member Functions

template<typename TElem >
void ass_dA_elem (LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
 computes the local defect and transfers it to the global discretization More...
 
template<typename TElem >
void ass_dM_elem (LocalVector &d, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
 computes the mass part of the defect of a time-dependent problem More...
 
template<size_t numEdges>
void ass_elem_stiffness (number perm, number cond, number S[2][numEdges][2][numEdges])
 composes the stiffness matrix of the stationary problem More...
 
template<typename TElem >
void ass_JA_elem (LocalMatrix &J, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
 transfers the precomputed local stiffness matrix to the global discretization More...
 
template<typename TElem >
void ass_JM_elem (LocalMatrix &J, const LocalVector &u, GridObject *elem, const position_type vCornerCoords[])
 computes the mass matrix of a time-dependent problem More...
 
template<typename TElem >
void ass_rhs_elem (LocalVector &d, GridObject *elem, const position_type vCornerCoords[])
 computes the right-hand side due to the generator currents More...
 
template<typename TElem >
void finish_element_loop ()
 finalizes the loop over the elements: clear the source More...
 
template<typename TElem >
void prepare_element (const LocalVector &u, GridObject *elem, ReferenceObjectID roid, const position_type vCornerCoords[])
 prepares a given element for assembling: computes the discretization of the rot-rot operator More...
 
template<typename TElem >
void prepare_element_loop (ReferenceObjectID roid, int si)
 prepares the loop over the elements: checks whether the parameters are set, ... More...
 
virtual void prepare_setting (const std::vector< LFEID > &vLfeID, bool bNonRegular)
 check type of the grid and the trial space More...
 
void register_all_loc_discr_funcs ()
 registers the local assembler functions for all the elements and dimensions More...
 
template<typename TElem >
void register_loc_discr_func ()
 registers the local assembler functions for a given element More...
 

Private Attributes

number m_conductivity
 the electric conductivity in the subdomain More...
 
number m_omega
 frequency \(\omega\) for the discretization More...
 
number m_permeability
 the magnetic permeability in the subdomain More...
 
tGeneratorCurrentm_pSsJG
 the source active in the current (assembled) subset More...
 
number m_rot_rot_M [maxNumEdges][maxNumEdges]
 local mass matrix of the rot-rot operator More...
 
number m_rot_rot_S [maxNumEdges][maxNumEdges]
 local stiffness matrix of the rot-rot operator More...
 
ConstSmartPtr< EMaterial< domain_type > > m_spSubsetData
 parameters of the materials in the domain More...
 
std::vector< tGeneratorCurrentm_vJG
 array of all the sources (generator currents) More...
 

Static Private Attributes

static const int dim = base_type::dim
 world dimension More...
 
static const size_t maxNumEdges = element_list_traits<typename domain_traits<dim>::DimElemList>::maxEdges
 max. number of the edges of the full-dimensional elements in the domain More...
 

Additional Inherited Members

- Public Types inherited from ug::IElemDisc< TDomain >
typedef IElemAssembleFuncs< IElemDisc< TDomain >, TDomain > assemble_base_type
 
typedef IElemError< TDomain > base_type
 
typedef TDomain domain_type
 
typedef IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain > estimator_base_type
 
- Public Types inherited from IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >
typedef TDomain domain_type
 
typedef IElemDisc< TDomain > leaf_type
 
- Public Types inherited from ug::IElemError< class >
typedef TDomain domain_type
 
typedef IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain > estimator_base_type
 
- Public Types inherited from ug::IElemDiscBase< class >
typedef TDomain domain_type
 
typedef TDomain::position_type position_type
 
- Public Types inherited from IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain >
typedef TDomain domain_type
 
typedef IElemDisc< TDomain > leaf_type
 
- Static Public Attributes inherited from ug::IElemDisc< TDomain >
static const int dim
 
- Static Public Attributes inherited from IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >
static const int dim
 
- Static Public Attributes inherited from ug::IElemError< class >
static const int dim
 
- Static Public Attributes inherited from ug::IElemDiscBase< class >
static const int dim
 
- Static Public Attributes inherited from IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain >
static const int dim
 
- Static Public Attributes inherited from ug::Electromagnetism::EddyCurrentTraits
static const size_t _Im_ = 1
 index of the imaginary part in the grid functions More...
 
static const size_t _Re_ = 0
 index of the real part in the grid functions More...
 
- Protected Member Functions inherited from ug::IElemDisc< TDomain >
void clear_add_fct ()
 
void clear_add_fct (ReferenceObjectID id)
 
void set_default_add_fct ()
 
- Protected Member Functions inherited from IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >
void clear_add_fct ()
 
void clear_add_fct (ReferenceObjectID id)
 
void remove_add_def_A_elem_fct (ReferenceObjectID id)
 
void remove_add_def_A_expl_elem_fct (ReferenceObjectID id)
 
void remove_add_def_M_elem_fct (ReferenceObjectID id)
 
void remove_add_jac_A_elem_fct (ReferenceObjectID id)
 
void remove_add_jac_M_elem_fct (ReferenceObjectID id)
 
void remove_add_rhs_elem_fct (ReferenceObjectID id)
 
void remove_fsh_elem_loop_fct (ReferenceObjectID id)
 
void remove_fsh_timestep_elem_fct (ReferenceObjectID id)
 
void remove_fsh_timestep_fct (size_t algebra_id)
 
void remove_prep_elem_fct (ReferenceObjectID id)
 
void remove_prep_elem_loop_fct (ReferenceObjectID id)
 
void remove_prep_timestep_elem_fct (ReferenceObjectID id)
 
void remove_prep_timestep_fct (size_t algebra_id)
 
void set_add_def_A_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_add_def_A_expl_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_add_def_M_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_add_jac_A_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_add_jac_M_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_add_rhs_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_default_add_fct ()
 
void set_fsh_elem_loop_fct (ReferenceObjectID id, TAssFunc func)
 
void set_fsh_timestep_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_fsh_timestep_fct (size_t algebra_id, TAssFunc func)
 
void set_prep_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_prep_elem_loop_fct (ReferenceObjectID id, TAssFunc func)
 
void set_prep_timestep_elem_fct (ReferenceObjectID id, TAssFunc func)
 
void set_prep_timestep_fct (size_t algebra_id, TAssFunc func)
 
- Protected Member Functions inherited from ug::IElemError< class >
void clear_add_fct ()
 
void clear_add_fct (ReferenceObjectID id)
 
- Protected Member Functions inherited from ug::IElemDiscBase< class >
virtual void approximation_space_changed ()
 
void set_function_pattern (ConstSmartPtr< FunctionPattern > fctPatt)
 
void update_function_index_mapping ()
 
- Protected Member Functions inherited from IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain >
void clear_add_fct ()
 
void clear_add_fct (ReferenceObjectID id)
 
void remove_compute_err_est_A_elem (ReferenceObjectID id)
 
void remove_compute_err_est_M_elem (ReferenceObjectID id)
 
void remove_compute_err_est_rhs_elem (ReferenceObjectID id)
 
void remove_fsh_err_est_elem_loop (ReferenceObjectID id)
 
void remove_prep_err_est_elem (ReferenceObjectID id)
 
void remove_prep_err_est_elem_loop (ReferenceObjectID id)
 
void set_compute_err_est_A_elem (ReferenceObjectID id, TAssFunc func)
 
void set_compute_err_est_M_elem (ReferenceObjectID id, TAssFunc func)
 
void set_compute_err_est_rhs_elem (ReferenceObjectID id, TAssFunc func)
 
void set_default_add_fct ()
 
void set_fsh_err_est_elem_loop (ReferenceObjectID id, TAssFunc func)
 
void set_prep_err_est_elem (ReferenceObjectID id, TAssFunc func)
 
void set_prep_err_est_elem_loop (ReferenceObjectID id, TAssFunc func)
 
- Protected Attributes inherited from ug::IElemDisc< TDomain >
std::vector< SmartPtr< IElemDiscModifier< TDomain > > > m_spElemModifier
 
- Protected Attributes inherited from IElemAssembleFuncs< IElemDisc< TDomain >, TDomain >
ReferenceObjectID m_roid
 
- Protected Attributes inherited from ug::IElemDiscBase< class >
bool m_bStationaryForced
 
FunctionGroup m_fctGrp
 
FunctionIndexMapping m_fctIndexMap
 
LocalVectorTimeSeriesm_pLocalVectorTimeSeries
 
SmartPtr< ApproximationSpace< TDomain > > m_spApproxSpace
 
ConstSmartPtr< FunctionPatternm_spFctPattern
 
size_t m_timePoint
 
std::vector< std::string > m_vFct
 
std::vector< IDataImport< dim > * > m_vIImport
 
std::vector< numberm_vScaleMass
 
std::vector< numberm_vScaleStiff
 
std::vector< std::string > m_vSubset
 
- Protected Attributes inherited from IElemEstimatorFuncs< IElemDisc< TDomain >, TDomain >
ReferenceObjectID m_roid
 
SmartPtr< IErrEstData< TDomain > > m_spErrEstData
 

Detailed Description

template<typename TDomain, typename TAlgebra>
class ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >

FE-discretization of the time-harmonic E-based formulation of the eddy current model.

This class implements the local FE-discretization of the time-harmonic E-based formulation of the eddy current model using the linear Nedelec-Type-1 (Whitney-1) elements. The model is formulated as follows:

\begin{eqnarray*} \mathbf{rot} \mu^{-1} \mathbf{rot} \mathbf{E} + i \omega \sigma \mathbf{E} = - i \omega \mathbf{J}_G, \end{eqnarray*}

where

  • \( \mathbf{E} \) (unknown, complex-valued scalar) electric field
  • \( \mu \) (given) magnetic permeability (a real number)
  • \( \sigma \) (given) conductivity (a real number)
  • \( \omega \) (given) frequency (a real number)
  • \( \mathbf{J}_G \) (given) generator current (non-zero only in the subdomain with \( \sigma = 0 \))
  • \( i = \sqrt{-1} \)
Remarks
Note that every Nedelec shape function requires a scalar (complex) dof value. These functions are thus vector-valued for the scalar dofs, each dof storing two doubles (one for the real and one for the imaginary part of the value).

References:

  • O. Sterz. Modellierung und Numerik zeitharmonischer Wirbelstromprobleme in 3D. PhD thesis, Heidelberg University, 2003.
  • A. Bossavit. Computational Electromagnetism. Academic Press (Boston), 1998
Template Parameters
TDomainDomain type
TAlgebraAlgebra type

Member Typedef Documentation

◆ base_type

template<typename TDomain , typename TAlgebra >
typedef IElemDisc<TDomain> ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::base_type
private

base class type

◆ domain_type

template<typename TDomain , typename TAlgebra >
typedef base_type::domain_type ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::domain_type
private

domain type

◆ position_type

template<typename TDomain , typename TAlgebra >
typedef base_type::position_type ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::position_type
private

position type

◆ TGridFunction

template<typename TDomain , typename TAlgebra >
typedef GridFunction<TDomain, TAlgebra> ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::TGridFunction
private

type of grid functions (used for the sources)

◆ this_type

template<typename TDomain , typename TAlgebra >
typedef EddyCurrent_E_Nedelec<TDomain,TAlgebra> ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::this_type
private

own type

Constructor & Destructor Documentation

◆ EddyCurrent_E_Nedelec()

template<typename TDomain , typename TAlgebra >
ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::EddyCurrent_E_Nedelec ( const char *  functions,
ConstSmartPtr< EMaterial< domain_type > >  spSubsetData,
number  frequency 
)

class constructor

Parameters
functionsgrid function names
spSubsetDatamaterial data object
frequencyfrequency \(\omega\)

References ug::IElemDiscBase< class >::num_fct(), ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::register_all_loc_discr_funcs(), and UG_THROW.

Member Function Documentation

◆ ass_dA_elem()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::ass_dA_elem ( LocalVector d,
const LocalVector u,
GridObject elem,
const position_type  vCornerCoords[] 
)
private

computes the local defect and transfers it to the global discretization

◆ ass_dM_elem()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::ass_dM_elem ( LocalVector d,
const LocalVector u,
GridObject elem,
const position_type  vCornerCoords[] 
)
private

computes the mass part of the defect of a time-dependent problem

This is a stationary problem, so there is no mass matrix in this sence. For this, the following function does nothing.

◆ ass_elem_stiffness()

template<typename TDomain , typename TAlgebra >
template<size_t numEdges>
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::ass_elem_stiffness ( number  perm,
number  cond,
number  S[2][numEdges][2][numEdges] 
)
private

composes the stiffness matrix of the stationary problem

composes the local stiffness matrix

This function composes the local stiffness matrix of the stationary problem by combining the stiffness and the mass matrix of the discretization of the rot-rot operator. Following O. Sterz, the system has the following block form:

\begin{eqnarray*} \left ( \begin{array} {cc} - S_h (\mu^{-1}) & \omega M^{(1)}_h (\sigma) \\ \omega M^{(1)}_h (\sigma) & S_h (\mu^{-1}) \end{array} \right ) \left ( \begin{array} {c} \mathbf{Re} \mathbf{E}_h \\ \mathbf{Im} \mathbf{E}_h \end{array} \right ) = \left ( \begin{array} {c} - \omega M^{(1)}_h (1) \mathbf{Im} \mathbf{J}_{Gh} \\ - \omega M^{(1)}_h (1) \mathbf{Re} \mathbf{J}_{Gh} \end{array} \right ), \end{eqnarray*}

where

  • \(S_h (\mu^{-1})\) the stiffness matrix of the \(\mathbf{rot}\mu^{-1}\mathbf{rot}\),
  • \( M^{(1)}_h (\sigma)\) the mass matrix of the \(\sigma I\) operator,
  • \(\mathbf{E}_h\) the unknown electric field grid function (represented using the Nedelec element),
  • \(\mathbf{J}_{Gh}\) the generator current (represented using the Nedelec element).
Remarks
For the solvability of the problem, the generator current must be numerically weakly divergence free (cf. the documentation for the class NedelecLoopCurrent for details).
Parameters
permthe magnetic permeability
condthe electric conductivity
Sfor the composed matrix

◆ ass_JA_elem()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::ass_JA_elem ( LocalMatrix J,
const LocalVector u,
GridObject elem,
const position_type  vCornerCoords[] 
)
private

transfers the precomputed local stiffness matrix to the global discretization

◆ ass_JM_elem()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::ass_JM_elem ( LocalMatrix J,
const LocalVector u,
GridObject elem,
const position_type  vCornerCoords[] 
)
private

computes the mass matrix of a time-dependent problem

This is a stationary problem, so there is no mass matrix in this sence. For this, the following function does nothing.

◆ ass_rhs_elem()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::ass_rhs_elem ( LocalVector b,
GridObject elem,
const position_type  vCornerCoords[] 
)
private

computes the right-hand side due to the generator currents

The right-hand side has the form

\begin{eqnarray*} \left ( \begin{array} {c} - \omega M^{(1)}_h \mathbf{Im} \mathbf{J}_{Gh} \\ - \omega M^{(1)}_h \mathbf{Re} \mathbf{J}_{Gh} \end{array} \right ), \end{eqnarray*}

where

  • \( M^{(1)}_h\) the mass matrix of the \(I\) (identity) operator,
  • \(\mathbf{J}_{Gh}\) the generator current (represented using the Nedelec element).
Remarks
For the solvability of the problem, the generator current must be numerically weakly divergence free.

References ug::DoFRef().

◆ finish_element_loop()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::finish_element_loop
private

finalizes the loop over the elements: clear the source

◆ prepare_element()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::prepare_element ( const LocalVector u,
GridObject elem,
ReferenceObjectID  roid,
const position_type  vCornerCoords[] 
)
private

prepares a given element for assembling: computes the discretization of the rot-rot operator

Parameters
ulocal solution at the dofs associated with elem
elemelement to prepare
vCornerCoordscoordinates of the corners of the element

References get(), and ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::local_stiffness_and_mass().

◆ prepare_element_loop()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::prepare_element_loop ( ReferenceObjectID  roid,
int  si 
)
private

prepares the loop over the elements: checks whether the parameters are set, ...

Parameters
roidonly elements with this roid are looped over
siand only in this subdomain

References UG_THROW.

◆ prepare_setting()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::prepare_setting ( const std::vector< LFEID > &  vLfeID,
bool  bNonRegular 
)
privatevirtual

check type of the grid and the trial space

Implements ug::IElemDiscBase< class >.

References ug::LFEID::NEDELEC, type(), and UG_THROW.

◆ register_all_loc_discr_funcs()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::register_all_loc_discr_funcs
private

registers the local assembler functions for all the elements and dimensions

Referenced by ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::EddyCurrent_E_Nedelec().

◆ register_loc_discr_func()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::register_loc_discr_func
private

registers the local assembler functions for a given element

Referenced by ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::RegisterLocalDiscr::operator()().

◆ set_generator_current() [1/2]

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::set_generator_current ( SmartPtr< TGridFunction spgfJG,
const char *  cmp 
)
inline

adds a generator current item \( \mathbf{J}_{G,h} \) to the discretization

Parameters
spgfJGgrid function with the data
cmpgrid function components

References ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::set_generator_current().

◆ set_generator_current() [2/2]

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::set_generator_current ( SmartPtr< TGridFunction spgfJG,
const char *  cmp,
const char *  ss_names = NULL 
)

adds a generator current item \( \mathbf{J}_{G,h} \) to the discretization

Adds a generator current item \( \mathbf{J}_{G,h} \) to the discretization.

Remarks
Note that for the solvability of the problem, the correct specification of the subsets is important.
See also
NedelecLoopCurrent
Parameters
spgfJGgrid function with the data
cmpgrid function components
ss_namesnames of the subsets where the current is defined (NULL for "everywhere")

References ug::SubsetGroup::add(), dim, SmartPtr< TGridFunction >::invalid(), ug::LFEID::NEDELEC, ug::RemoveWhitespaceFromString(), ug::TokenizeString(), UG_CATCH_THROW, and UG_THROW.

Referenced by ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::set_generator_current().

Member Data Documentation

◆ dim

template<typename TDomain , typename TAlgebra >
const int ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::dim = base_type::dim
staticprivate

world dimension

◆ m_conductivity

template<typename TDomain , typename TAlgebra >
number ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::m_conductivity
private

the electric conductivity in the subdomain

◆ m_omega

template<typename TDomain , typename TAlgebra >
number ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::m_omega
private

frequency \(\omega\) for the discretization

◆ m_permeability

template<typename TDomain , typename TAlgebra >
number ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::m_permeability
private

the magnetic permeability in the subdomain

◆ m_pSsJG

template<typename TDomain , typename TAlgebra >
tGeneratorCurrent* ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::m_pSsJG
private

the source active in the current (assembled) subset

◆ m_rot_rot_M

template<typename TDomain , typename TAlgebra >
number ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::m_rot_rot_M[maxNumEdges][maxNumEdges]
private

local mass matrix of the rot-rot operator

◆ m_rot_rot_S

template<typename TDomain , typename TAlgebra >
number ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::m_rot_rot_S[maxNumEdges][maxNumEdges]
private

local stiffness matrix of the rot-rot operator

◆ m_spSubsetData

template<typename TDomain , typename TAlgebra >
ConstSmartPtr<EMaterial<domain_type> > ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::m_spSubsetData
private

parameters of the materials in the domain

◆ m_vJG

template<typename TDomain , typename TAlgebra >
std::vector<tGeneratorCurrent> ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::m_vJG
private

array of all the sources (generator currents)

◆ maxNumEdges

template<typename TDomain , typename TAlgebra >
const size_t ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::maxNumEdges = element_list_traits<typename domain_traits<dim>::DimElemList>::maxEdges
staticprivate

max. number of the edges of the full-dimensional elements in the domain


The documentation for this class was generated from the following files: