Plugins
|
FE-discretization of the time-harmonic E-based formulation of the eddy current model. More...
#include <eddy_current_e_nedelec.h>
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 | |
void | set_generator_current (SmartPtr< TGridFunction > spgfJG, const char *cmp) |
adds a generator current item \( \mathbf{J}_{G,h} \) to the discretization | |
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 | |
![]() | |
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) | |
![]() | |
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[]) |
TLeaf & | 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 () |
![]() | |
IElemError (const char *functions, const char *subsets) | |
IElemError (const std::vector< std::string > &vFct, const std::vector< std::string > &vSubset) | |
![]() | |
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 FunctionGroup & | function_group () const |
ConstSmartPtr< FunctionPattern > | function_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 LocalVectorTimeSeries * | local_time_solutions () const |
const FunctionIndexMapping & | map () 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_type & | subset_handler () |
const TDomain::subset_handler_type & | subset_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 () |
![]() | |
TLeaf & | 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 | |
typedef base_type::domain_type | domain_type |
domain type | |
typedef base_type::position_type | position_type |
position type | |
typedef GridFunction< TDomain, TAlgebra > | TGridFunction |
type of grid functions (used for the sources) | |
typedef EddyCurrent_E_Nedelec< TDomain, TAlgebra > | this_type |
own type | |
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 | |
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 | |
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 | |
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 | |
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 | |
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 | |
template<typename TElem > | |
void | finish_element_loop () |
finalizes the loop over the elements: clear the source | |
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 | |
template<typename TElem > | |
void | prepare_element_loop (ReferenceObjectID roid, int si) |
prepares the loop over the elements: checks whether the parameters are set, ... | |
virtual void | prepare_setting (const std::vector< LFEID > &vLfeID, bool bNonRegular) |
check type of the grid and the trial space | |
void | register_all_loc_discr_funcs () |
registers the local assembler functions for all the elements and dimensions | |
template<typename TElem > | |
void | register_loc_discr_func () |
registers the local assembler functions for a given element | |
Private Attributes | |
number | m_conductivity |
the electric conductivity in the subdomain | |
number | m_omega |
frequency \(\omega\) for the discretization | |
number | m_permeability |
the magnetic permeability in the subdomain | |
tGeneratorCurrent * | m_pSsJG |
the source active in the current (assembled) subset | |
number | m_rot_rot_M [maxNumEdges][maxNumEdges] |
local mass matrix of the rot-rot operator | |
number | m_rot_rot_S [maxNumEdges][maxNumEdges] |
local stiffness matrix of the rot-rot operator | |
ConstSmartPtr< EMaterial< domain_type > > | m_spSubsetData |
parameters of the materials in the domain | |
std::vector< tGeneratorCurrent > | m_vJG |
array of all the sources (generator currents) | |
Static Private Attributes | |
static const int | dim = base_type::dim |
world dimension | |
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 | |
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
References:
TDomain | Domain type |
TAlgebra | Algebra type |
|
private |
base class type
|
private |
domain type
|
private |
position type
|
private |
type of grid functions (used for the sources)
|
private |
own type
ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::EddyCurrent_E_Nedelec | ( | const char * | functions, |
ConstSmartPtr< EMaterial< domain_type > > | spSubsetData, | ||
number | frequency | ||
) |
class constructor
functions | grid function names |
spSubsetData | material data object |
frequency | frequency \(\omega\) |
References ug::IElemDiscBase< typename TDomain >::num_fct(), ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::register_all_loc_discr_funcs(), and UG_THROW.
|
private |
computes the local defect and transfers it to the global discretization
|
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.
|
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
perm | the magnetic permeability |
cond | the electric conductivity |
S | for the composed matrix |
|
private |
transfers the precomputed local stiffness matrix to the global discretization
|
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.
|
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
References ug::DoFRef().
|
private |
finalizes the loop over the elements: clear the source
|
private |
prepares a given element for assembling: computes the discretization of the rot-rot operator
u | local solution at the dofs associated with elem |
elem | element to prepare |
vCornerCoords | coordinates of the corners of the element |
References ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::local_stiffness_and_mass().
|
private |
prepares the loop over the elements: checks whether the parameters are set, ...
roid | only elements with this roid are looped over |
si | and only in this subdomain |
References UG_THROW.
|
privatevirtual |
check type of the grid and the trial space
Implements ug::IElemDiscBase< typename TDomain >.
References ug::LFEID::NEDELEC, type(), and UG_THROW.
|
private |
registers the local assembler functions for all the elements and dimensions
Referenced by ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::EddyCurrent_E_Nedelec().
|
private |
registers the local assembler functions for a given element
Referenced by ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::RegisterLocalDiscr::operator()().
|
inline |
adds a generator current item \( \mathbf{J}_{G,h} \) to the discretization
spgfJG | grid function with the data |
cmp | grid function components |
References ug::Electromagnetism::EddyCurrent_E_Nedelec< TDomain, TAlgebra >::set_generator_current().
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.
spgfJG | grid function with the data |
cmp | grid function components |
ss_names | names of the subsets where the current is defined (NULL for "everywhere") |
References ug::SubsetGroup::add(), SmartPtr< typename T, template< class TT > class FreePolicy >::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().
|
staticprivate |
world dimension
|
private |
the electric conductivity in the subdomain
|
private |
frequency \(\omega\) for the discretization
|
private |
the magnetic permeability in the subdomain
|
private |
the source active in the current (assembled) subset
|
private |
local mass matrix of the rot-rot operator
|
private |
local stiffness matrix of the rot-rot operator
|
private |
parameters of the materials in the domain
|
private |
array of all the sources (generator currents)
|
staticprivate |
max. number of the edges of the full-dimensional elements in the domain