Plugins
ug::Electromagnetism Namespace Reference

Classes

class  CalcVolPowerElemHelperClass
 Helper class for the integration of the power. More...
 
class  CalcVolPowerHelperClass
 Helper class for the computation of the power of the electromagnetic field. More...
 
class  ComPol_MemCopy
 
class  ComPol_MemOp
 
class  ComputeElemFluxHelper
 Helper class for computation of the flux in elements. (Helper for ComputeFlux.) More...
 
class  ComputeElemFluxHelper< TGridFunc, RegularEdge >
 Helper class for computation of the flux in elements. (Helper for ComputeFlux.) More...
 
class  ComputeFluxHelper
 Helper class for the loop over all the elements in the computation of the flux. (Helper for ComputeFlux.) More...
 
class  EddyCurrent_E_Nedelec
 FE-discretization of the time-harmonic E-based formulation of the eddy current model. More...
 
class  EddyCurrentCurlEDependentCmpUserData
 
class  EddyCurrentHeat
 
class  EddyCurrentImBofEUserData
 
class  EddyCurrentReBofEUserData
 
class  EddyCurrentTraits
 Auxiliary class defining some important constants. More...
 
class  EMaterial
 Class for subdomain-dependent data for the E-based formulated problems. More...
 
class  EMDirichlet
 Common interface to get the Dirichlet boundary conditions. More...
 
struct  Functionality
 
class  NedelecCurlData
 
class  NedelecDirichletBC
 Dirichlet BC for a rot-rot operator. More...
 
class  NedelecGridFunctionData
 
class  NedelecInterpolation
 Interpolation of the Nedelec dofs and their curls. More...
 
class  NedelecInterpolation< TDomain, 2, 2 >
 A specialization of NedelecInterpolation for 2d,. More...
 
class  NedelecInterpolation< TDomain, 3, 3 >
 A specialization of NedelecInterpolation for 3d,. More...
 
class  NedelecLoopCurrent
 Class for computation of loop currents. More...
 
class  NedelecProject
 Projection procedure class. More...
 
class  NedelecProlongationMatrixHelper
 
class  NedelecProlongationMatrixHelper< TDomain, TAlgebra, RegularEdge >
 
class  NedelecSigmaEData
 
class  NedelecT1_LDisc
 Tool kit for the Whitney-1 (Nedelec) based FE discretization of the rot-rot operators. More...
 
class  NedelecT1_LDisc< TDomain, Tetrahedron >
 Specialization of NedelecT1_LDisc for tetrahedra. More...
 
class  NedelecT1_LDisc< TDomain, Triangle >
 Specialization of NedelecT1_LDisc for triangles. More...
 
class  NedelecT1_LDisc_forSimplex
 Helper class for the specialization of NedelecT1_LDisc for simplices (triangles and tetrahedrons) More...
 
class  NedelecTransfer
 Class of the prolongation and the restriction of the Nedelec DoFs. More...
 
class  TimeHarmonicNedelecHybridSmoother
 The hybrid smoother by R. Hiptmair. More...
 

Functions

template<typename TGridFunc >
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. More...
 
template<typename TGridFunc >
number calc_magnetic_flux (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 flux[2])
 Computation of the magnetic flux through windings of a coil. More...
 
template<typename TGridFunc >
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) More...
 
template<typename TGridFunc >
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 a coil. More...
 
template<typename TGridFunc >
void CalcMagneticFlux (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 magnetic flux through windings of a coil. More...
 
template<typename TGridFunc >
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. More...
 
template<typename TGridFunc >
void ComputeFlux (SmartPtr< TGridFunc > spGF, const char *fct_names, const char *vol_subsets, const char *face_subsets)
 prints the flux through a given surface in the shell More...
 
template<typename TGridFunc >
number ComputeFlux (TGridFunc *pGF, size_t fct, SubsetGroup &volSSG, SubsetGroup &faceSSG)
 returns the flux through a given surface More...
 
template<typename TDomain >
number ComputeNedelecDoF (UserData< MathVector< TDomain::dim >, TDomain::dim > *function, Edge *edge, const typename TDomain::position_accessor_type &aaPos, const int si, number time)
 
template<typename TGridFunction >
void ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, size_t fct, const SubsetGroup &ssGrp, number time)
 
template<int dim>
void get_ave_vector (const size_t n_vec, const MathVector< dim > *vec, MathVector< dim > &ave)
 
template<typename TArray , typename TLayout , typename TOp >
void MemAllReduce (TArray *pMem, const TLayout &masterLayout, const TLayout &slaveLayout, pcl::InterfaceCommunicator< TLayout > *pCom=NULL)
 
template<typename TPotGridFunc , typename TGridFunc >
void NedelecGradPotential (SmartPtr< TPotGridFunc > spPotGF, SmartPtr< TGridFunc > spGF, const char *cmp)
 
template<typename TPotGridFunc , typename TGridFunc >
void NedelecGradPotential (SmartPtr< TPotGridFunc > spPotGF, SmartPtr< TGridFunc > spGF, size_t fct)
 
template<typename TPotGridFunc >
void SetSubsetVertVal (SmartPtr< TPotGridFunc > spPotGF, const char *subsets, number value)
 
template<typename TPotGridFunc >
void SetSubsetVertVal (SmartPtr< TPotGridFunc > spPotGF, const SubsetGroup &ssGrp, number value)
 
template<typename TGridFunction >
void ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, const char *cmp, const char *subsets, number time)
 
template<typename TGridFunction >
void ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, const char *cmp, number time)
 
template<typename TGridFunction >
void ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, const char *cmp, const char *subsets)
 
template<typename TGridFunction >
void ComputeNedelecDoFs (SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spFunction, SmartPtr< TGridFunction > spGridFct, const char *cmp)
 

Function Documentation

◆ calc_EMF()

template<typename TGridFunc >
number ug::Electromagnetism::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.

This function computes the magnetic flux through windings of a cylindric coil and returns the electomotive force due to the varying magnetic flux. The flux is a complex number (as well as the electric field that is used as the base of the computation). The function returns the total area of all the windings.

Parameters
[in]gfEgrid function with the electric field
[in]fctfunction components (Re and Im)
[in]ss_grpsubsets of the kerner of the coil
[in]Normalnormal to the planes of the windings
[in]base_pntpoint on the plane of the 1st winding
[in]n_pntnumber of the windings
[in]d_pntthickness of the windings
[out]emfthe computed electromotive force (Re and Im)

References ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim >::curl_flux(), dim, ug::DoFRef(), ug::LFEID::NEDELEC, normal, ug::SubsetGroup::size(), UG_THROW, and ug::VecLength().

Referenced by CalcEMF().

◆ calc_magnetic_flux()

template<typename TGridFunc >
number ug::Electromagnetism::calc_magnetic_flux ( 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  flux[2] 
)

Computation of the magnetic flux through windings of a coil.

Parameters
[in]gfEgrid function with the electric field
[in]fctfunction components (Re and Im)
[in]ss_grpsubsets of the kerner of the coil
[in]Normalnormal to the planes of the windings
[in]base_pntpoint on the plane of the 1st winding
[in]n_pntnumber of the windings
[in]d_pntthickness of the windings
[out]fluxthe computed flux (Re and Im)

◆ calc_power()

template<typename TGridFunc >
void ug::Electromagnetism::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)

Parameters
pJGGFgrid function of the Nedelec-DoFs of the generator current \(\mathbf{J}_G\)
JG_fctindices of the Re and Im parts in the grid function for \(\mathbf{J}_G\)
JG_ssg(full-dim.) subsets where \(\mathbf{J}_G\) is defined (non-zero and in the kernel)
pEGFgrid function of the Nedelec-DoFs of the electric field \(\mathbf{E}\)
E_fctindices of the Re and Im parts in the grid function for \(\mathbf{E}\)
powto add the integral

References dim.

Referenced by CalcPower().

◆ CalcEMF()

template<typename TGridFunc >
void ug::Electromagnetism::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 a coil.

Parameters
[in]spGFgrid function with the electric field
[in]cmpsnames of the components of the grid function (for Re and Im)
[in]subsetssubsets where to compute the flux
[in]Normalnormal to the planes of the windings
[in]base_pntpoint on the plane of the 1st winding
[in]n_pntnumber of the windings
[in]d_pntthickness of the windings

References ug::SubsetGroup::add(), calc_EMF(), dim, SmartPtr< class, FreePolicy >::get(), ug::RemoveWhitespaceFromString(), ug::TokenizeString(), UG_CATCH_THROW, UG_LOG, and UG_THROW.

◆ CalcMagneticFlux()

template<typename TGridFunc >
void ug::Electromagnetism::CalcMagneticFlux ( 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 magnetic flux through windings of a coil.

Parameters
[in]spGFgrid function with the electric field
[in]cmpsnames of the components of the grid function (for Re and Im)
[in]subsetssubsets where to compute the flux
[in]Normalnormal to the planes of the windings
[in]base_pntpoint on the plane of the 1st winding
[in]n_pntnumber of the windings
[in]d_pntthickness of the windings

◆ CalcPower()

template<typename TGridFunc >
void ug::Electromagnetism::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.

Parameters
[in]spJGGFgrid function with the generator current
[in]JG_cmpsnames of the components of the grid function (for Re and Im)
JG_ss(full-dim.) subsets where \(\mathbf{J}_G\) is defined (non-zero and in the kernel)
[in]spEGFgrid function with the electric field
[in]E_cmpsnames of the components of the grid function (for Re and Im)

References ug::SubsetGroup::add(), calc_power(), SmartPtr< class, FreePolicy >::get(), ug::RemoveWhitespaceFromString(), ug::SubsetGroup::set_subset_handler(), ug::TokenizeString(), UG_CATCH_THROW, UG_LOG, and UG_THROW.

◆ ComputeFlux() [1/2]

template<typename TGridFunc >
void ug::Electromagnetism::ComputeFlux ( SmartPtr< TGridFunc >  spGF,
const char *  fct_name,
const char *  vol_subsets,
const char *  face_subsets 
)

prints the flux through a given surface in the shell

Computation of the flux of a vector field (given by the Whitney-1 elements) through a low-dimensional subset (surface). The function prints the value of the flux in the shell.

Parameters
spGFgrid function of the Nedelec-DoFs of the vector field
fct_namename of the function in the grid function
vol_subsetsfull-dim. subsets (adjacent to the surface) to indicate the negative direction
face_subsetsthe surface (the low-dim. subsets)

References ug::SubsetGroup::add(), ComputeFlux(), SmartPtr< class, FreePolicy >::get(), ug::RemoveWhitespaceFromString(), ug::SubsetGroup::set_subset_handler(), ug::FunctionGroup::size(), ug::TokenizeString(), UG_CATCH_THROW, UG_LOG, and UG_THROW.

◆ ComputeFlux() [2/2]

template<typename TGridFunc >
number ug::Electromagnetism::ComputeFlux ( TGridFunc *  pGF,
size_t  fct,
SubsetGroup volSSG,
SubsetGroup faceSSG 
)

returns the flux through a given surface

Computation of the flux of a vector field (given by the Whitney-1 elements) through a low-dimensional subset (surface). The function returns the value of the flux.

Parameters
pGFgrid function of the Nedelec-DoFs of the vector field
fctindex of the function in the grid function
volSSGfull-dim. subsets (adjacent to the surface) to indicate the negative direction
faceSSGthe surface (the low-dim. subsets)

References dim, ug::LFEID::NEDELEC, and UG_THROW.

Referenced by ComputeFlux().

◆ ComputeNedelecDoF()

template<typename TDomain >
number ug::Electromagnetism::ComputeNedelecDoF ( UserData< MathVector< TDomain::dim >, TDomain::dim > *  function,
Edge edge,
const typename TDomain::position_accessor_type aaPos,
const int  si,
number  time 
)
inline

Computes the value for the DoF associated with a given edge. Uses the trapezoid quadrature rule.

Parameters
functionthe function to encode
edgethe edge to compute the DoF for
aaPoscoordinates of the vertices
sisubset index
timethe time argument

References ug::EdgeVertices::vertex().

◆ ComputeNedelecDoFs() [1/5]

template<typename TGridFunction >
void ug::Electromagnetism::ComputeNedelecDoFs ( SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > >  spFunction,
SmartPtr< TGridFunction >  spGridFct,
const char *  cmp 
)

Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function. This function parses the string specifications of the function components and the subsets.

Parameters
spFunctionthe function to encode
spGridFctthe grid function to store the result
cmpname of the component

References ComputeNedelecDoFs().

◆ ComputeNedelecDoFs() [2/5]

template<typename TGridFunction >
void ug::Electromagnetism::ComputeNedelecDoFs ( SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > >  spFunction,
SmartPtr< TGridFunction >  spGridFct,
const char *  cmp,
const char *  subsets 
)

Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function. This function parses the string specifications of the function components and the subsets.

Parameters
spFunctionthe function to encode
spGridFctthe grid function to store the result
cmpname of the component
subsetsnames of the subsets (may be NULL meaning "all")

References ComputeNedelecDoFs().

◆ ComputeNedelecDoFs() [3/5]

template<typename TGridFunction >
void ug::Electromagnetism::ComputeNedelecDoFs ( SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > >  spFunction,
SmartPtr< TGridFunction >  spGridFct,
const char *  cmp,
const char *  subsets,
number  time 
)

Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function. This function parses the string specifications of the function components and the subsets.

Parameters
spFunctionthe function to encode
spGridFctthe grid function to store the result
cmpname of the component
subsetsnames of the subsets (may be NULL meaning "all")
timethe time argument

References ug::SubsetGroup::add(), ug::SubsetGroup::add_all(), ComputeNedelecDoFs(), ug::LFEID::NEDELEC, ug::TokenizeString(), and UG_THROW.

◆ ComputeNedelecDoFs() [4/5]

template<typename TGridFunction >
void ug::Electromagnetism::ComputeNedelecDoFs ( SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > >  spFunction,
SmartPtr< TGridFunction >  spGridFct,
const char *  cmp,
number  time 
)

Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function. This function parses the string specifications of the function components and the subsets.

Parameters
spFunctionthe function to encode
spGridFctthe grid function to store the result
cmpname of the component
timethe time argument

References ComputeNedelecDoFs().

◆ ComputeNedelecDoFs() [5/5]

template<typename TGridFunction >
void ug::Electromagnetism::ComputeNedelecDoFs ( SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > >  spFunction,
SmartPtr< TGridFunction >  spGridFct,
size_t  fct,
const SubsetGroup ssGrp,
number  time 
)

Computes the values of the DoFs of the Whitney-1 (Nedelec-type-1) basis functions for a given vector function.

Parameters
spFunctionthe function to encode
spGridFctthe grid function to store the result
fctindex of the function component in the grid function
ssGrpsubset group where to encode
timethe time argument (for the function to encode)

References ug::DoFRef(), ug::SubsetGroup::size(), UG_ASSERT, UG_LOG, and UG_THROW.

Referenced by ComputeNedelecDoFs().

◆ get_ave_vector()

template<int dim>
void ug::Electromagnetism::get_ave_vector ( const size_t  n_vec,
const MathVector< dim > *  vec,
MathVector< dim > &  ave 
)
inline

Computes the arithm. average of given MathVectors.

Parameters
[in]n_vecnumber of the vectors in the array
[in]vecarray of the vectors
[out]avethe average

◆ MemAllReduce()

template<typename TArray , typename TLayout , typename TOp >
void ug::Electromagnetism::MemAllReduce ( TArray *  pMem,
const TLayout &  masterLayout,
const TLayout &  slaveLayout,
pcl::InterfaceCommunicator< TLayout > *  pCom = NULL 
)

Syncronizes a memory array using a given reduction operation.

Template Parameters
TArraytype of the array
TLayouttype of the parallel layout
TOpclass of the reduction operation
Parameters
pMemthe array to syncronize
masterLayoutthe master layout
slaveLayoutthe slave layout
pCom(optionally) the communicator

References pcl::InterfaceCommunicator< class >::communicate(), pcl::InterfaceCommunicator< class >::receive_data(), and pcl::InterfaceCommunicator< class >::send_data().

◆ NedelecGradPotential() [1/2]

template<typename TPotGridFunc , typename TGridFunc >
void ug::Electromagnetism::NedelecGradPotential ( SmartPtr< TPotGridFunc >  spPotGF,
SmartPtr< TGridFunc >  spGF,
const char *  cmp 
)

Computes the Whitney representation of the gradient of a vertex-centered scalar grid function. (This is the version with the symbolic function name.)

Parameters
spPotGFthe potential to compute the gradient for
spGFto save the Whitney form of the gradient
cmpname of the component

References NedelecGradPotential(), and UG_THROW.

◆ NedelecGradPotential() [2/2]

template<typename TPotGridFunc , typename TGridFunc >
void ug::Electromagnetism::NedelecGradPotential ( SmartPtr< TPotGridFunc >  spPotGF,
SmartPtr< TGridFunc >  spGF,
size_t  fct 
)

Computes the Whitney representation of the gradient of a vertex-centered scalar grid function.

Parameters
spPotGFthe potential to compute the gradient for
spGFto save the Whitney form of the gradient
fctindex of the function in spGF

References ug::DoFRef(), ug::LFEID::NEDELEC, UG_THROW, and ug::EdgeVertices::vertex().

Referenced by NedelecGradPotential().

◆ SetSubsetVertVal() [1/2]

template<typename TPotGridFunc >
void ug::Electromagnetism::SetSubsetVertVal ( SmartPtr< TPotGridFunc >  spPotGF,
const char *  subsets,
number  value 
)

Sets a scalar potential field to a given value at vertices of a given subset, and to 0 for all other subsets.

Parameters
spPotGFthe potential to compute the gradient for
subsetsnames of the subsets (may be NULL meaning "all")
valuevalue to set in the subset group

References ug::SubsetGroup::add(), ug::SubsetGroup::add_all(), SetSubsetVertVal(), and ug::TokenizeString().

◆ SetSubsetVertVal() [2/2]

template<typename TPotGridFunc >
void ug::Electromagnetism::SetSubsetVertVal ( SmartPtr< TPotGridFunc >  spPotGF,
const SubsetGroup ssGrp,
number  value 
)

Sets a scalar potential field to a given value at vertices of a given subset, and to 0 for all other subsets.

Parameters
spPotGFthe potential to compute the gradient for
ssGrpsubset group where to set
valuevalue to set in the subset group

References ug::SubsetGroup::size(), and UG_THROW.

Referenced by SetSubsetVertVal().