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

Class for computation of loop currents. More...

#include <nedelec_source.h>

Classes

class  AuxLaplaceLocAss
 Class for local assembling of the auxiliary Laplace operator. More...
 
struct  GetFluxOfPotential
 Helper class for computation of the flux of potential over the cut. More...
 
class  LocLaplaceA
 Helper class for assembling the local stiffness matrix of Laplacian in various routines. More...
 
struct  MarkSourceEdges
 Helper class for marking the edges in the source. More...
 
class  OutOfSource
 Constraint that sets the problem to 0-identity out of the source. More...
 
struct  TSrcData
 Structure for keeping electric current data. More...
 
class  ZeroAverage
 Postprocessor to eliminate the kernel parts in the solution. More...
 

Public Types

typedef TAlgebra algebra_type
 Type of algebra (for the Nedelec-element-based grid functions) More...
 
typedef TDomain domain_type
 Type of Domain. More...
 
typedef TDomain::grid_type grid_type
 Type of Grid: More...
 
typedef TAlgebra::matrix_type matrix_type
 Type of Vector (for the Nedelec-element-based grid functions) More...
 
typedef TDomain::position_type position_type
 position type More...
 
typedef GridFunction< TDomain, TPotAlgebrapot_gf_type
 Grid function type for the potential space. More...
 
typedef TPotAlgebra::matrix_type pot_matrix_type
 Matrix type for the potential space. More...
 
typedef TPotAlgebra::vector_type pot_vector_type
 Vector type for the potential space. More...
 
typedef domain_type::subset_handler_type subset_handler_type
 Type of subset handler. More...
 
typedef CPUAlgebra TPotAlgebra
 The auxiliary algebra type for the space of the potential. (Note: It should be scalar.) More...
 
typedef TAlgebra::vector_type vector_type
 Type of Vector (for the Nedelec-element-based grid functions) More...
 

Public Member Functions

void compute (SmartPtr< GridFunction< TDomain, TAlgebra > > sp_u)
 Computation of the source (only if 'init'ialized): distributes the potential. More...
 
 NedelecLoopCurrent (const char *ssNames, const char *posSsNames, const char *cutSsNames, SmartPtr< ApproximationSpace< TDomain > > vertApproxSpace, SmartPtr< ILinearOperatorInverse< pot_vector_type > > potSolver)
 Constructor. More...
 
void set (const char *fctNames, number I)
 Setting the electric current value. More...
 
std::string subsets ()
 Returns the subsets where the source is defined. More...
 
SmartPtr< ITransferPostProcess< TDomain, TPotAlgebra > > zero_average ()
 Returns the postprocess routine to eliminate the kernel part of the potential. More...
 

Static Public Attributes

static const int WDim = TDomain::dim
 world dimention More...
 

Private Types

typedef AChar a_edge_flag_type
 Type of the attachment and its accessor for flags. More...
 
typedef ABool a_vert_flag_type
 
typedef Grid::EdgeAttachmentAccessor< a_edge_flag_typeaa_edge_flag_type
 
typedef Grid::VertexAttachmentAccessor< a_vert_flag_typeaa_vert_flag_type
 
typedef NedelecLoopCurrent< TDomain, TAlgebra > this_type
 This type. More...
 

Private Member Functions

void compute_potential (pot_gf_type &pot_u)
 Computes the potential of the source. More...
 
void distribute_source_potential (const DoFDistribution &vertDD, pot_vector_type &src_pot, DoFDistribution &edgeDD, size_t func, number value, vector_type &src_field)
 Computes the gradient of the potential. More...
 
template<typename TElem >
void get_flux_of_pot (const TDomain &domain, const pot_vector_type &pot, const DoFDistribution &vertDD, number &flux)
 Computes the flux of of the gradient of the potential over the cut (for one type of elements) More...
 
template<typename TElem >
void mark_source_edges (const DoFDistribution &edgeDD, aa_edge_flag_type &in_source)
 Marks edges that belong to the loop source domain (for one type of elements) More...
 

Private Attributes

SubsetGroup m_allSsGrp
 Group of all the subsets of the source. More...
 
std::string m_allSsNames
 Names of all the subsets of the source. More...
 
SmartPtr< DomainDiscretization< TDomain, TPotAlgebra > > m_auxLaplaceAss
 Global discretization of the auxiliary equations. More...
 
SmartPtr< AssembledLinearOperator< TPotAlgebra > > m_auxLaplaceOp
 Matrix of the discretization. More...
 
SmartPtr< AuxLaplaceLocAssm_auxLocLaplace
 Local discretization of the auxiliary equations. More...
 
SubsetGroup m_cutSsGrp
 Group of the surfaces of the cut of the loop. More...
 
std::string m_cutSsNames
 Names of the surfaces of the cut of the loop. More...
 
SmartPtr< OutOfSourcem_outOfSource
 Extension of the matrices and vectors to the whole domain. More...
 
SubsetGroup m_posSsGrp
 Group of the subsets of the positive direction. More...
 
std::string m_posSsNames
 Names of the subsets of the positive direction. More...
 
SmartPtr< ILinearOperatorInverse< pot_vector_type > > m_potSolver
 Solver for the auxliliary equations. More...
 
SmartPtr< ApproximationSpace< TDomain > > m_spVertApproxSpace
 Vertex-centered approximation space. More...
 
std::vector< TSrcDatam_vSrcData
 Array of the electric current data. More...
 
SmartPtr< ZeroAveragem_zeroAverage
 Elimination of the kernal parts in the solution. More...
 

Detailed Description

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

Class for computation of loop currents.

Class for computation of weakly divergence-free loop currents that are non-zero only in a given subdomain \(T\) (i.e. in a coil).

The range (image) of the stiffness matrix of the discrete \(\mathbf{rot} \, \mathbf{rot}\) operator is the subspace of all functions othrogonal to all the gradients \( \mathbf{G} \psi \), where \( \mathbf{G} \) is the incidence matrix between the vertices and the edges, and \( \psi \) is any vertex-centered grid function. Thus, for the solvability of the discretized system (for ex. of the E-based formulation of the eddy-current model), its right-hand side \(M^{(1)} J \) ( \(M^{(1)}\) being the mass matrix of the Whitney-1 shapes) must satisfy

\begin{eqnarray*} ( \mathbf{G}^T M^{(1)} J, \psi ) = ( M^{(1)} J, \mathbf{G} \psi ) = 0 \end{eqnarray*}

This is really true only in insulators. In conductors, a mass matrix is added to the \(\mathbf{rot} \, \mathbf{rot}\) operator, so that the kernel is trivial. We assume that a coil is modelled as a zero-conductivity medium with the source, so that \(T\) is completely contained in insulators. Thus the former condition implies

\begin{eqnarray*} \mathbf{G}_{i}^T M^{(1)} J = 0, \end{eqnarray*}

(with \( J = 0 \) out of \( T \)) i.e. the weak divergence of the rhs must be zero.

This class implements a computation of such a weakly divergence-free current \(J\) in a given set of insulators \(T\) constituting (topologically) a ring (torus) by specifying a cut with the jump of the potential and a part adjacent to this cut that specifies the positive direction of the current. The current is represented as \(J = \hat{\mathbf{G}} \phi\) and the potential \(\phi\) is computed as the solution of a discretized Laplace equation with the special boundary conditions. Here, \(\hat{\mathbf{G}}\) is the matrix obtained from \(\mathbf{G}\) by replacing the rows

The computed \(J\) is rescaled so that the current over the cut is equal to the given value. Different values of the currents can be used for different function names: cf. the 'set' functions of this class.

Note that the mass matrix \(M^{(1)}\) is not the same as the mass matrix at the off-diagonal positions of the stiffness matrix of the time-harmonic problem. The matrix \(M^{(1)}\) mentioned here for the computation of the right-hand side is assembled only of the local mass matrices for the grid elements from \(T\). The resulting discretized system for \(\phi\) consists the discrete Laplace operator on \(T\), the discrete Neumann-0 (natural) boundary condition on \(\partial{T}\) and the identity matrix for all other DoFs. Note furthermore that the same \(M^{(1)}\) should be used in the assembling of the right-hand side of the discretized Maxwell equations.

References:

  • O. Sterz. Modellierung und Numerik zeitharmonischer Wirbelstromprobleme in 3D. PhD thesis, 2003.
Template Parameters
TDomainDomain type
TAlgebraAlgebra type

Member Typedef Documentation

◆ a_edge_flag_type

template<typename TDomain , typename TAlgebra >
typedef AChar ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::a_edge_flag_type
private

Type of the attachment and its accessor for flags.

◆ a_vert_flag_type

template<typename TDomain , typename TAlgebra >
typedef ABool ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::a_vert_flag_type
private

◆ aa_edge_flag_type

template<typename TDomain , typename TAlgebra >
typedef Grid::EdgeAttachmentAccessor<a_edge_flag_type> ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::aa_edge_flag_type
private

◆ aa_vert_flag_type

template<typename TDomain , typename TAlgebra >
typedef Grid::VertexAttachmentAccessor<a_vert_flag_type> ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::aa_vert_flag_type
private

◆ algebra_type

template<typename TDomain , typename TAlgebra >
typedef TAlgebra ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::algebra_type

Type of algebra (for the Nedelec-element-based grid functions)

◆ domain_type

template<typename TDomain , typename TAlgebra >
typedef TDomain ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::domain_type

Type of Domain.

◆ grid_type

template<typename TDomain , typename TAlgebra >
typedef TDomain::grid_type ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::grid_type

Type of Grid:

◆ matrix_type

template<typename TDomain , typename TAlgebra >
typedef TAlgebra::matrix_type ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::matrix_type

Type of Vector (for the Nedelec-element-based grid functions)

◆ position_type

template<typename TDomain , typename TAlgebra >
typedef TDomain::position_type ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::position_type

position type

◆ pot_gf_type

template<typename TDomain , typename TAlgebra >
typedef GridFunction<TDomain, TPotAlgebra> ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::pot_gf_type

Grid function type for the potential space.

◆ pot_matrix_type

template<typename TDomain , typename TAlgebra >
typedef TPotAlgebra::matrix_type ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::pot_matrix_type

Matrix type for the potential space.

◆ pot_vector_type

template<typename TDomain , typename TAlgebra >
typedef TPotAlgebra::vector_type ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::pot_vector_type

Vector type for the potential space.

◆ subset_handler_type

template<typename TDomain , typename TAlgebra >
typedef domain_type::subset_handler_type ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::subset_handler_type

Type of subset handler.

◆ this_type

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

This type.

◆ TPotAlgebra

template<typename TDomain , typename TAlgebra >
typedef CPUAlgebra ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::TPotAlgebra

The auxiliary algebra type for the space of the potential. (Note: It should be scalar.)

◆ vector_type

template<typename TDomain , typename TAlgebra >
typedef TAlgebra::vector_type ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::vector_type

Type of Vector (for the Nedelec-element-based grid functions)

Constructor & Destructor Documentation

◆ NedelecLoopCurrent()

template<typename TDomain , typename TAlgebra >
ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::NedelecLoopCurrent ( const char *  ssNames,
const char *  posSsNames,
const char *  cutSsNames,
SmartPtr< ApproximationSpace< TDomain > >  vertApproxSpace,
SmartPtr< ILinearOperatorInverse< pot_vector_type > >  potSolver 
)

Member Function Documentation

◆ compute()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::compute ( SmartPtr< GridFunction< TDomain, TAlgebra > >  sp_u)

◆ compute_potential()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::compute_potential ( pot_gf_type pot_u)
private

◆ distribute_source_potential()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::distribute_source_potential ( const DoFDistribution vertDD,
pot_vector_type src_pot,
DoFDistribution edgeDD,
size_t  func,
number  value,
vector_type src_field 
)
private

Computes the gradient of the potential.

Computes the gradient of the potential for one function.

Remark: The potential is considered to be properly scaled. Otherwise, the scaling should be included into the value of the source

Parameters
[in]vertDDthe vertex DD
[in]src_potpotential to distribute
[in]edgeDDthe edge DD
[in]funcindex of the function
[in]valuevalue of the source
[out]src_fieldthe computed source field

References ug::DoFDistribution::begin(), ug::DoFRef(), ug::DoFDistribution::end(), ug::DoFDistribution::inner_algebra_indices(), ug::DoFDistribution::inner_dof_indices(), ug::DoFDistribution::multi_grid(), PCL_RO_BOR, ug::SetAttachmentValues(), UG_THROW, and ug::EdgeVertices::vertex().

◆ get_flux_of_pot()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::get_flux_of_pot ( const TDomain &  domain,
const pot_vector_type pot,
const DoFDistribution vertDD,
number flux 
)
private

Computes the flux of of the gradient of the potential over the cut (for one type of elements)

Computes the flux of the gradient of the potential over the cut (for one type of elements)

Parameters
[in]domainthe domain
[in]potthe potential field
[in]vertDDthe vertex DD
[out]fluxthe flux to update

References ug::DoFDistribution::inner_algebra_indices(), ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_allSsGrp, ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_cutSsGrp, ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_posSsGrp, ug::SubsetGroup::remove(), ug::SubsetGroup::size(), ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::LocLaplaceA< TElem >::stiffness(), and UG_THROW.

◆ mark_source_edges()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::mark_source_edges ( const DoFDistribution edgeDD,
aa_edge_flag_type in_source 
)
private

Marks edges that belong to the loop source domain (for one type of elements)

Marks edges that belong to the loop source domain (for one type of elements). For edges of elements belonging to the source domain, the entries of 'in_source' are set to non-zero. (Other entries are not set.) For those entries that correspond to edges belonging to elements in the 'positive direction' subdomain and having one of the ends at the cut surface, the 1st and the 2nd bits of the entry is set to 1 (depending whether the beginning or the end lyies at the cut). For all other edges in the source, these bits are set to zero (so that only the 0th bit is 1).

Parameters
[in]edgeDDthe edge DD
[out]in_sourcethe flags to update

References ug::Grid::associated_elements(), ug::ISubsetHandler::get_subset_index(), ug::ISubsetHandler::grid(), ug::PointerConstArray< class >::size(), ug::DoFDistributionInfoProvider::subset_handler(), UG_ASSERT, and ug::EdgeVertices::vertex().

◆ set()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::set ( const char *  fctNames,
number  I 
)

Setting the electric current value.

Setting the electric current value

Parameters
fctNamesnames of the components
Ielectric current for the functions

◆ subsets()

template<typename TDomain , typename TAlgebra >
std::string ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::subsets ( )
inline

Returns the subsets where the source is defined.

References ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_allSsNames.

◆ zero_average()

template<typename TDomain , typename TAlgebra >
SmartPtr<ITransferPostProcess<TDomain, TPotAlgebra> > ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::zero_average ( )
inline

Returns the postprocess routine to eliminate the kernel part of the potential.

References ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_zeroAverage.

Member Data Documentation

◆ m_allSsGrp

template<typename TDomain , typename TAlgebra >
SubsetGroup ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_allSsGrp
private

◆ m_allSsNames

template<typename TDomain , typename TAlgebra >
std::string ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_allSsNames
private

◆ m_auxLaplaceAss

template<typename TDomain , typename TAlgebra >
SmartPtr<DomainDiscretization<TDomain, TPotAlgebra> > ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_auxLaplaceAss
private

Global discretization of the auxiliary equations.

This discretization object uses m_auxLocLaplace and m_outOfSource to assemble the Laplace equation and the jump right-hand side in the subsets of the source, as well as the identity equations in the other parts of the domain.

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

◆ m_auxLaplaceOp

template<typename TDomain , typename TAlgebra >
SmartPtr<AssembledLinearOperator<TPotAlgebra> > ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_auxLaplaceOp
private

Matrix of the discretization.

◆ m_auxLocLaplace

template<typename TDomain , typename TAlgebra >
SmartPtr<AuxLaplaceLocAss> ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_auxLocLaplace
private

Local discretization of the auxiliary equations.

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

◆ m_cutSsGrp

template<typename TDomain , typename TAlgebra >
SubsetGroup ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_cutSsGrp
private

◆ m_cutSsNames

template<typename TDomain , typename TAlgebra >
std::string ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_cutSsNames
private

Names of the surfaces of the cut of the loop.

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

◆ m_outOfSource

template<typename TDomain , typename TAlgebra >
SmartPtr<OutOfSource> ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_outOfSource
private

Extension of the matrices and vectors to the whole domain.

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

◆ m_posSsGrp

◆ m_posSsNames

template<typename TDomain , typename TAlgebra >
std::string ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_posSsNames
private

Names of the subsets of the positive direction.

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

◆ m_potSolver

template<typename TDomain , typename TAlgebra >
SmartPtr<ILinearOperatorInverse<pot_vector_type> > ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_potSolver
private

◆ m_spVertApproxSpace

template<typename TDomain , typename TAlgebra >
SmartPtr<ApproximationSpace<TDomain> > ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_spVertApproxSpace
private

◆ m_vSrcData

template<typename TDomain , typename TAlgebra >
std::vector<TSrcData> ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_vSrcData
private

Array of the electric current data.

◆ m_zeroAverage

template<typename TDomain , typename TAlgebra >
SmartPtr<ZeroAverage> ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::m_zeroAverage
private

Elimination of the kernal parts in the solution.

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

◆ WDim

template<typename TDomain , typename TAlgebra >
const int ug::Electromagnetism::NedelecLoopCurrent< TDomain, TAlgebra >::WDim = TDomain::dim
static

world dimention


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