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

The hybrid smoother by R. Hiptmair. More...

#include <hiptmair_hybrid_smoother.h>

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

Classes

struct  t_red_op_and
 "and" reduction operation class for the conductivity condition More...
 
struct  t_red_op_or
 "or" reduction operation class for the conductivity condition More...
 
struct  tEdgeInfo
 Structure of the storage for the information about the edge-vertex interconnections. More...
 

Public Types

typedef MatrixOperator< matrix_type, vector_typematrix_operator_type
 Matrix Operator type. More...
 
typedef TAlgebra::matrix_type matrix_type
 Matrix type. More...
 
typedef MatrixOperator< pot_matrix_type, pot_vector_typepot_matrix_operator_type
 Matrix Operator 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 GridFunction< TDomain, TAlgebra > TGridFunc
 Grid function type for the solution. More...
 
typedef CPUAlgebra TPotAlgebra
 The auxiliary algebra type for the space of the potential. (Note: It should be scalar.) More...
 
typedef GridFunction< TDomain, TPotAlgebraTPotGridFunc
 Grid function type for the potential. More...
 
typedef TAlgebra::vector_type vector_type
 Vector type. More...
 
- Public Types inherited from ug::ILinearIterator< class, class >
typedef Y codomain_function_type
 
typedef X domain_function_type
 
- Public Types inherited from ug::DebugWritingObject< TAlgebra >
typedef TAlgebra algebra_type
 
typedef TAlgebra::matrix_type matrix_type
 
typedef TAlgebra::vector_type vector_type
 
- Public Types inherited from VectorDebugWritingObject< TAlgebra::vector_type >
typedef TAlgebra::vector_type vector_type
 

Public Member Functions

bool apply (vector_type &c, const vector_type &d)
 Computes the correction. More...
 
bool apply_update_defect (vector_type &c, vector_type &d)
 Computes the correction and updates the defect d := d - L*c. More...
 
SmartPtr< ILinearIterator< vector_type, vector_type > > clone ()
 Clone the smoother by copying the data. More...
 
bool init (SmartPtr< ILinearOperator< vector_type > > J, const vector_type &u)
 Initialization using a matrix and a GridFunction (not merely a vector!) More...
 
bool init (SmartPtr< ILinearOperator< vector_type > > L)
 We cannot initialize without the geometry (this version of init cannot not work) More...
 
const char * name () const
 Returns the name. More...
 
void set_Dirichlet (SmartPtr< EMDirichlet< TDomain, TAlgebra > > spDirichlet)
 Sets the Dirichlet boundary. More...
 
void set_skip_edge_smoother (bool skip_edge)
 Skip flag the edge smoother. More...
 
void set_skip_vertex_smoother (bool skip_vertex)
 Skip flag the vertex smoother. More...
 
bool supports_parallel () const
 Currently returns false because the computation of the potentials is not purely parallel up to now. More...
 
 TimeHarmonicNedelecHybridSmoother (SmartPtr< ApproximationSpace< TDomain > > vertApproxSpace, SmartPtr< ILinearIterator< vector_type > > edgeSmoother, SmartPtr< ILinearIterator< pot_vector_type > > vertSmoother)
 Constructor setting the approx. spaces and the default subsmoothers. More...
 
 ~TimeHarmonicNedelecHybridSmoother ()
 Destructor. More...
 
- Public Member Functions inherited from ug::ILinearIterator< class, class >
virtual bool apply (Y &c, const X &d)=0
 
virtual bool apply_update_defect (Y &c, X &d)=0
 
virtual std::string config_string () const
 
SmartPtr< IDamping< X, Y > > damping ()
 
 ILinearIterator ()
 
 ILinearIterator (const ILinearIterator< X, Y > &parent)
 
virtual bool init (SmartPtr< ILinearOperator< Y, X > > J, const Y &u)=0
 
virtual bool init (SmartPtr< ILinearOperator< Y, X > > L)=0
 
void set_damp (number factor)
 
void set_damp (SmartPtr< IDamping< X, Y > > spScaling)
 
virtual ~ILinearIterator ()
 
- Public Member Functions inherited from ug::DebugWritingObject< TAlgebra >
SmartPtr< IDebugWriter< algebra_type > > debug_writer ()
 
ConstSmartPtr< IDebugWriter< algebra_type > > debug_writer () const
 
bool debug_writer_valid () const
 
 DebugWritingObject ()
 
 DebugWritingObject (const DebugWritingObject< algebra_type > &parent)
 
 DebugWritingObject (SmartPtr< IDebugWriter< algebra_type > > spDebugWriter)
 
virtual void set_debug (SmartPtr< IDebugWriter< algebra_type > > spDebugWriter)
 
virtual ~DebugWritingObject ()
 
- Public Member Functions inherited from VectorDebugWritingObject< TAlgebra::vector_type >
virtual void set_debug (SmartPtr< IVectorDebugWriter< vector_type > > spDebugWriter)
 
SmartPtr< IVectorDebugWriter< vector_type > > vector_debug_writer ()
 
ConstSmartPtr< IVectorDebugWriter< vector_type > > vector_debug_writer () const
 
bool vector_debug_writer_valid () const
 
 VectorDebugWritingObject ()
 
 VectorDebugWritingObject (SmartPtr< IVectorDebugWriter< vector_type > > spDebugWriter)
 
void write_debug (const vector_type &vec, const char *filename)
 
virtual ~VectorDebugWritingObject ()
 

Private Types

typedef TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra > this_type
 

Private Member Functions

void collect_edge_defect (const vector_type &d, pot_vector_type &potDefRe, pot_vector_type &potDefIm)
 Computes the vertex-centered defect \( d_{pot} = G^T d \). More...
 
void compute_GtMG ()
 Computes the product \( G^T M^{(1)}_h G \). More...
 
void compute_potential_matrix (const DoFDistribution *pEdgeDD, const DoFDistribution *pVertDD)
 Computes the matrix for the smoother in the potential space and marks the "conductive nodes". More...
 
void distribute_vertex_correction (pot_vector_type &potCorRe, pot_vector_type &potCorIm, vector_type &c)
 Adds the vertex-centered correction to the edge-centered one: More...
 
void get_edge_vert_correspondence (const DoFDistribution *pEdgeDD, const DoFDistribution *pVertDD)
 Gets the correspondence between the edges and the vertices. More...
 
bool zero_mass_entry (double val)
 Measure for numerically zero entries in the mass matrix. More...
 

Private Attributes

bool m_bInit
 Whether initialized. More...
 
bool m_bSkipEdge
 Needed mainly for debugging: Whether to skip one of the stages. More...
 
bool m_bSkipVertex
 
GridLevel m_GridLevel
 Grid level (or surface) where to smooth. More...
 
TPotGridFuncm_pPotCorIm
 Vertex-centered grid function for the Im-part of potential corrections (this is a GridFunction to allow geometry-dependent smoothers) More...
 
TPotGridFuncm_pPotCorRe
 Vertex-centered grid function for the Re-part of potential corrections (this is a GridFunction to allow geometry-dependent smoothers) More...
 
SmartPtr< EMDirichlet< TDomain, TAlgebra > > m_spDirichlet
 Dirichlet boundary. More...
 
SmartPtr< DoFDistributionm_spEdgeDD
 DoF distribution for the Nedelec elements. More...
 
SmartPtr< ILinearIterator< vector_type > > m_spEdgeSmoother
 Smoother for the edge dofs. More...
 
SmartPtr< pot_matrix_operator_typem_spPotMat
 Vertex-centered matrix for the potential. More...
 
SmartPtr< matrix_operator_typem_spSysMat
 Edge-centered matrix of the original equation. More...
 
SmartPtr< ApproximationSpace< TDomain > > m_spVertApproxSpace
 Approximation space for the potential (vertex-centered) space. More...
 
SmartPtr< DoFDistributionm_spVertDD
 Level DoF distribution for the vertex centered grid func. More...
 
SmartPtr< ILinearIterator< pot_vector_type > > m_spVertSmoother
 Smoother for the vertex dofs. More...
 
VariableArray1< tEdgeInfom_vEdgeInfo
 Storage for the information about the edge-vertex interconnections. More...
 

Additional Inherited Members

- Protected Member Functions inherited from ug::DebugWritingObject< TAlgebra >
void enter_debug_writer_section (const char *secDir)
 
void enter_debug_writer_section (std::string secDir)
 
void leave_debug_writer_section ()
 
void write_debug (const matrix_type &mat, const char *filename)
 
void write_debug (const matrix_type &mat, std::string name)
 
- Protected Member Functions inherited from VectorDebugWritingObject< TAlgebra::vector_type >
void enter_vector_debug_writer_section (const char *secDir)
 
void enter_vector_debug_writer_section (std::string secDir)
 
void leave_vector_debug_writer_section ()
 
void print_debugger_message (const char *msg)
 
void print_debugger_message (std::string msg)
 
virtual void write_debug (const vector_type &vec, std::string name)
 
- Protected Attributes inherited from ug::ILinearIterator< class, class >
SmartPtr< IDamping< X, Y > > m_spDamping
 
- Protected Attributes inherited from ug::DebugWritingObject< TAlgebra >
SmartPtr< IDebugWriter< algebra_type > > m_spDebugWriter
 
- Protected Attributes inherited from VectorDebugWritingObject< TAlgebra::vector_type >
SmartPtr< IVectorDebugWriter< vector_type > > m_spVectorDebugWriter
 

Detailed Description

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

The hybrid smoother by R. Hiptmair.

This class implements the hybrid smoother proposed by R. Hiptmair for the discretizations based on the Nedelec element. The smoother is based on the approximate Helmholtz decomponstion of the operator. This implementation is constructed for the discretizations of the time-harmonic equations. The class computes the matrix for the discrete potential and calls two further smoothers specified by the user: one in the space of the edge elements (Whitney-1 or Nedelec forms) and one in the vertex-centered potential space (i.e. the Whitney-0 or Lagrange-1 elements). Accordingly, this smoother requires 2 approximation spaces: one for the edge dofs and one for the vertex dofs.

References:

  • R. Hiptmair. Multigrid method for Maxwell’s equations. SIAM J. Numer. Anal., 36(1): pp. 204–225 (1998), DOI 10.1137/S0036142997326203.
  • R. Beck, P. Deuflhard, R. Hiptmair, R. Hoppe und B. Wohlmuth. Adaptive multilevel methods for edge element discretizations of Maxwell’s equations. Surveys on Mathematics for Industry, 8(3-4): pp. 271–312 (1999).
Template Parameters
TDomainType of Domain
TAlgebraType of Algebra

Member Typedef Documentation

◆ matrix_operator_type

template<typename TDomain , typename TAlgebra >
typedef MatrixOperator<matrix_type, vector_type> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::matrix_operator_type

Matrix Operator type.

◆ matrix_type

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

Matrix type.

◆ pot_matrix_operator_type

template<typename TDomain , typename TAlgebra >
typedef MatrixOperator<pot_matrix_type, pot_vector_type> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::pot_matrix_operator_type

Matrix Operator type for the potential space.

◆ pot_matrix_type

template<typename TDomain , typename TAlgebra >
typedef TPotAlgebra::matrix_type ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< 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::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::pot_vector_type

Vector type for the potential space.

◆ TGridFunc

template<typename TDomain , typename TAlgebra >
typedef GridFunction<TDomain, TAlgebra> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::TGridFunc

Grid function type for the solution.

◆ this_type

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

◆ TPotAlgebra

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

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

◆ TPotGridFunc

template<typename TDomain , typename TAlgebra >
typedef GridFunction<TDomain, TPotAlgebra> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::TPotGridFunc

Grid function type for the potential.

◆ vector_type

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

Vector type.

Constructor & Destructor Documentation

◆ TimeHarmonicNedelecHybridSmoother()

template<typename TDomain , typename TAlgebra >
ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::TimeHarmonicNedelecHybridSmoother ( SmartPtr< ApproximationSpace< TDomain > >  vertApproxSpace,
SmartPtr< ILinearIterator< vector_type > >  edgeSmoother,
SmartPtr< ILinearIterator< pot_vector_type > >  vertSmoother 
)
inline

◆ ~TimeHarmonicNedelecHybridSmoother()

Member Function Documentation

◆ apply()

template<typename TDomain , typename TAlgebra >
bool ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::apply ( vector_type c,
const vector_type d 
)

Computes the correction.

This method applies the iterator operator, i.e. invertes the approximate matrix. The defect d remains unchanged.

The idea is to apply first the edge-centered smoother (specified as the argument of the constructor of the class). To the correction obtained in this way, the vertex-centered correction is added. The latter is obtained by the application of the vertex-centered smoother to the system

\begin{eqnarray*} \left ( \begin{array} {cc} 0 & A_{pot} \\ A_{pot} & 0 \end{array} \right ) \left ( \begin{array} {c} \mathbf{Re} c_{pot} \\ \mathbf{Im} c_{pot} \end{array} \right ) = \left ( \begin{array} {c} - \mathbf{Re} d_{pot} \\ \mathbf{Im} d_{pot} \end{array} \right ), \end{eqnarray*}

where \( A_{pot} = G^T M^{(1)}_h G \) and \( d_{pot} = G^T d \). To add the vertex-centered correction, obtained in this way, to the edge-centered one, one sets: \( c := c_{edge} + G c_{pot} \).

Parameters
[in]ddefect
[out]ccorrection
Returns
bool success flag

References ug::PST_ADDITIVE, ug::PST_CONSISTENT, ParallelVector< Vector< double > >::set_layouts(), ParallelVector< Vector< double > >::set_storage_type(), UG_LOG, and UG_THROW.

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

◆ apply_update_defect()

◆ clone()

◆ collect_edge_defect()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::collect_edge_defect ( const vector_type d,
pot_vector_type potDefRe,
pot_vector_type potDefIm 
)
private

Computes the vertex-centered defect \( d_{pot} = G^T d \).

Computes the defect in the potential space: \( d_{pot} = G^T d \), where \( d \) is the defect of the original system. This is done only for the 'conductive' vertices. At the other vertices, \( d_{pot} \) is set to 0.

Parameters
doriginal (edge-centered) defect
potDefRereal part of \( d_{pot} \)
potDefImimaginary part of \( d_{pot} \)

References ug::BlockRef(), ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::conductive_vrt_0(), ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::conductive_vrt_1(), ParallelVector< Vector< double > >::set(), and ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::vrt_index.

◆ compute_GtMG()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::compute_GtMG
private

Computes the product \( G^T M^{(1)}_h G \).

Performs the multiplication \( G^T M^{(1)}_h G \). Furthermore, marks the 'conductive' vertices, i.e. the vertices that lie on the closure of the conductive subset (i.e. fills the m_vConductiveVertex array). 'Conductive' vertices are ends of the edges whose entries of the mass matrix \( M^{(1)}_h\) are non-zero. Any other vertex, which is not an end of such an edge, is considered as 'non-conductive'.

For the non-conductive nodes, the function sets the row of the potential stiffness matrix to the identity.

This function uses the m_vEdgeInfo-array to get the correspondence between the edges.

References ug::VariableArray1< class >::at(), ug::BlockRef(), ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::is_Dirichlet(), ug::VariableArray1< class >::resize(), ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::set_conductive_vrt_0(), ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::set_conductive_vrt_1(), ug::SetDirichletRow(), UG_THROW, and ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::vrt_index.

◆ compute_potential_matrix()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::compute_potential_matrix ( const DoFDistribution pEdgeDD,
const DoFDistribution pVertDD 
)
private

Computes the matrix for the smoother in the potential space and marks the "conductive nodes".

Computes the matrix \(A_{pot}\) for the smoother in the potential space and marks the "conductive nodes" for the computation of the correction.

\begin{eqnarray*} A_{pot} = G^T M^{(1)}_h G, \end{eqnarray*}

where \(G\) is the note-to-edge gradient incidence matrix and \(M^{(1)}_h\) is the mass matrix of the Whitney-1-forms, taken from the (2, 1)-position of the block-form of the system matrix of the time-harmonic system:

\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 ) \end{eqnarray*}

(this is the edge-centered matrix).

Template Parameters
TDomaintype of domain in the class
TAlgebratype of algebra in the class
Parameters
pEdgeDDedge-centered DoF distribution of the grid functions
pVertDDvertex-centered DoF distribution of the grid functions

References ug::DoFDistribution::layouts(), ug::DoFDistribution::num_indices(), UG_CATCH_THROW, and UG_THROW.

◆ distribute_vertex_correction()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::distribute_vertex_correction ( pot_vector_type potCorRe,
pot_vector_type potCorIm,
vector_type c 
)
private

Adds the vertex-centered correction to the edge-centered one:

Updates the edge-centered correction: \( c := c_{edge} + G c_{pot} \).

Parameters
potCorRereal part of the potential correction \( c_{pot} \)
potCorImimaginary part of the potential correction \( c_{pot} \)
cfinal (edge-centered) correction

References ug::BlockRef(), and ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::vrt_index.

◆ get_edge_vert_correspondence()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::get_edge_vert_correspondence ( const DoFDistribution pEdgeDD,
const DoFDistribution pVertDD 
)
private

Gets the correspondence between the edges and the vertices.

Gets the correspondence between the edge DoFs and the vertex DoFs (i.e. fills the tEdgeInfo-structures in m_vEdgeInfo). Furthermore, this function marks edges at Dirichlet boundaries.

Template Parameters
TDomaintype of domain in the class
TAlgebratype of algebra in the class
Parameters
pEdgeDDedge-centered DoF distribution of the grid functions
pVertDDvertex-centered DoF distribution of the grid functions

References ug::DoFDistribution::begin(), ug::DoFDistribution::end(), ug::DoFDistribution::inner_algebra_indices(), ug::DoFDistribution::num_indices(), ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::set_init(), ug::SubsetGroup::size(), ug::DoFDistributionInfoProvider::subset_handler(), UG_THROW, ug::EdgeVertices::vertex(), and ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::tEdgeInfo::vrt_index.

◆ init() [1/2]

template<typename TDomain , typename TAlgebra >
bool ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::init ( SmartPtr< ILinearOperator< vector_type > >  J,
const vector_type u 
)

Initialization using a matrix and a GridFunction (not merely a vector!)

This function initializes the smoother with the matrix J and the geometry read from the GridFunction u. Note that u should be a GridFunction, not merely a vector.

Parameters
[in]Jthe matrix of the time-harmonic Nedelec discretization
[in]uthe grid function containing the grid info (the values are not used)
Returns
bool success flag

References ug::PST_CONSISTENT, and UG_THROW.

◆ init() [2/2]

template<typename TDomain , typename TAlgebra >
bool ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::init ( SmartPtr< ILinearOperator< vector_type > >  L)
inline

We cannot initialize without the geometry (this version of init cannot not work)

References ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::name(), and UG_THROW.

◆ name()

◆ set_Dirichlet()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::set_Dirichlet ( SmartPtr< EMDirichlet< TDomain, TAlgebra > >  spDirichlet)
inline

Sets the Dirichlet boundary.

Parameters
spDirichletthe Dirichlet BC object

References ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_spDirichlet.

◆ set_skip_edge_smoother()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::set_skip_edge_smoother ( bool  skip_edge)
inline

◆ set_skip_vertex_smoother()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::set_skip_vertex_smoother ( bool  skip_vertex)
inline

◆ supports_parallel()

template<typename TDomain , typename TAlgebra >
bool ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::supports_parallel ( ) const
inlinevirtual

Currently returns false because the computation of the potentials is not purely parallel up to now.

Implements ug::ILinearIterator< class, class >.

◆ zero_mass_entry()

template<typename TDomain , typename TAlgebra >
bool ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::zero_mass_entry ( double  val)
inlineprivate

Measure for numerically zero entries in the mass matrix.

Member Data Documentation

◆ m_bInit

template<typename TDomain , typename TAlgebra >
bool ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_bInit
private

Whether initialized.

◆ m_bSkipEdge

template<typename TDomain , typename TAlgebra >
bool ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_bSkipEdge
private

◆ m_bSkipVertex

◆ m_GridLevel

template<typename TDomain , typename TAlgebra >
GridLevel ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_GridLevel
private

Grid level (or surface) where to smooth.

◆ m_pPotCorIm

template<typename TDomain , typename TAlgebra >
TPotGridFunc* ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_pPotCorIm
private

Vertex-centered grid function for the Im-part of potential corrections (this is a GridFunction to allow geometry-dependent smoothers)

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

◆ m_pPotCorRe

template<typename TDomain , typename TAlgebra >
TPotGridFunc* ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_pPotCorRe
private

Vertex-centered grid function for the Re-part of potential corrections (this is a GridFunction to allow geometry-dependent smoothers)

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

◆ m_spDirichlet

◆ m_spEdgeDD

template<typename TDomain , typename TAlgebra >
SmartPtr<DoFDistribution> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_spEdgeDD
private

DoF distribution for the Nedelec elements.

◆ m_spEdgeSmoother

◆ m_spPotMat

template<typename TDomain , typename TAlgebra >
SmartPtr<pot_matrix_operator_type> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_spPotMat
private

Vertex-centered matrix for the potential.

◆ m_spSysMat

template<typename TDomain , typename TAlgebra >
SmartPtr<matrix_operator_type> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_spSysMat
private

◆ m_spVertApproxSpace

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

◆ m_spVertDD

template<typename TDomain , typename TAlgebra >
SmartPtr<DoFDistribution> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_spVertDD
private

Level DoF distribution for the vertex centered grid func.

◆ m_spVertSmoother

◆ m_vEdgeInfo

template<typename TDomain , typename TAlgebra >
VariableArray1<tEdgeInfo> ug::Electromagnetism::TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra >::m_vEdgeInfo
private

Storage for the information about the edge-vertex interconnections.


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