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

Projection procedure class. More...

#include <nedelec_project.h>

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

Classes

class  AuxLaplaceLocAss
 Class for local assembling of the auxiliary Laplace operator. More...
 
class  AuxLaplaceRHS
 Constraint that assembles the rhs und the bc for the auxiliary problems. More...
 
struct  ClearDivInConductors
 Helper class for assembling the weak divergence operator. More...
 
struct  IntegrateDivDVF
 Helper class for computation of the charges of the DVFs. More...
 
class  LocLaplaceA
 Helper class for assembling the local stiffness matrix of Laplacian in various routines. More...
 
struct  MarkCondVert
 Helper class for setting the base conductor indices to vertices. More...
 
struct  WeakDiv
 Helper class for the computation the weak divergence. 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 Types inherited from ug::IPProcessVector< TAlgebra::vector_type >
typedef TVector vector_type
 

Public Member Functions

void apply (GridFunction< TDomain, TAlgebra > &u, const char *fct_names)
 Performs the projection of given functions. More...
 
void apply (GridFunction< TDomain, TAlgebra > &u, const FunctionGroup &fctGrp)
 Performs the projection of given functions. More...
 
void apply (vector_type &vec)
 Performs the projection of all functions in a grid function. More...
 
void compute_div (SmartPtr< GridFunction< TDomain, TAlgebra > > sp_u, const char *u_fct_name, SmartPtr< GridFunction< TDomain, TPotAlgebra > > sp_div)
 Computes the weak divergence in insulators. More...
 
 NedelecProject (SmartPtr< EMaterial< TDomain > > emMatherial, SmartPtr< ApproximationSpace< TDomain > > vertApproxSpace, SmartPtr< ILinearOperatorInverse< pot_vector_type > > potSolver)
 Constructor. More...
 
void set_dampDVFs (bool val)
 Sets whether to damp Dirichlet vector fields. More...
 
void set_Dirichlet (SmartPtr< EMDirichlet< TDomain, TAlgebra > > spDirichlet)
 Sets the Dirichlet boundary. More...
 
 ~NedelecProject ()
 Destructor. More...
 
- Public Member Functions inherited from ug::IPProcessVector< TAlgebra::vector_type >
virtual void apply (vector_type &v)=0
 
virtual ~IPProcessVector ()
 

Static Public Attributes

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

Private Types

typedef AInt a_vert_cond_type
 Type of the attachment and its accessor the base conductor indices. More...
 
typedef Grid::VertexAttachmentAccessor< a_vert_cond_typeaa_vert_cond_type
 
typedef NedelecProject< TDomain, TAlgebra > this_type
 This type. More...
 

Private Member Functions

void alloc_DVFs (const TDomain &domain, pot_gf_type &aux_rhs)
 Allocates memory for the DVFs associated with the ungrounded conductors. More...
 
void assemble_div (const TDomain &domain, const DoFDistribution &edgeDD, const vector_type &u, size_t fct, const DoFDistribution &vertDD, pot_vector_type &div, DenseVector< VariableArray1< number > > &charge)
 Assembles the weak divergence operator on the whole grid. More...
 
template<typename TElem >
void clear_div_in_conductors (const DoFDistribution &vertDD, pot_vector_type &div, DenseVector< VariableArray1< number > > &charge)
 Sets the div operator to 0 in conductors (for all elements of the same type) More...
 
void compute_DVF_potential_coeffs (const TDomain &domain, DoFDistribution &vertDD)
 Computes the potential coefficients. More...
 
void compute_DVFs (pot_gf_type &aux_rhs)
 Computes the Dirichlet vector fields (DVFs) More...
 
void damp_DVFs (pot_vector_type &cor, const DenseVector< VariableArray1< number > > &charge)
 Damps the Dirichlet vector fields (DVFs) More...
 
void distribute_cor (const TDomain &domain, const DoFDistribution &edgeDD, vector_type &u, size_t fct, const DoFDistribution &vertDD, const pot_vector_type &cor)
 Updates the grid function by the potential correction. More...
 
template<typename TElem >
void integrate_div_DVF_elem_type (const TDomain &domain, const DoFDistribution &vertDD, const aa_vert_cond_type &vert_base_cond, DenseMatrix< VariableArray2< number > > &C)
 Integrates div E over boundaries of conductors for elements of the same type. More...
 
template<typename TElem >
void mark_cond_vert_elem_type (DoFDistribution &vertDD, aa_vert_cond_type &vert_base_cond)
 Sets the base conductor indices for every vertex (for all elements of the same type) More...
 
void project_func (const SmartPtr< TDomain > &domain, const SmartPtr< DoFDistribution > &edgeDD, vector_type &u, size_t fct, const SmartPtr< DoFDistribution > &vertDD, pot_gf_type &aux_rhs, pot_gf_type &aux_cor)
 Projects one function (i.e. performs the main action): More...
 
template<typename TElem >
void weak_div_elem_type (const TDomain &domain, const DoFDistribution &edgeDD, const vector_type &u, size_t fct, const DoFDistribution &vertDD, pot_vector_type &div)
 Computes the weak divergence (for all elements of the same type) More...
 

Private Attributes

SmartPtr< DomainDiscretization< TDomain, TPotAlgebra > > m_auxLaplaceAss
 Global discretization of the auxiliary equations. More...
 
SmartPtr< AssembledLinearOperator< TPotAlgebra > > m_auxLaplaceOp
 Matrix of the discretization. More...
 
SmartPtr< AuxLaplaceRHSm_auxLaplaceRHS
 Rhs and BC for the auxiliary problem. More...
 
SmartPtr< AuxLaplaceLocAssm_auxLocLaplace
 Local discretization of the auxiliary equations. More...
 
bool m_bDampDVFs
 Whether to damp the Dirichlet vector fields (DVFs) More...
 
std::vector< pot_gf_type * > m_DVF_phi
 Dirichlet vector fields (DVFs) More...
 
DenseMatrix< VariableArray2< number > > m_potCoeff
 Potential coefficients: More...
 
SmartPtr< ILinearOperatorInverse< pot_vector_type > > m_potSolver
 Solver for the auxliliary equations. More...
 
SmartPtr< EMDirichlet< TDomain, TAlgebra > > m_spDirichlet
 Dirichlet boundary. More...
 
SmartPtr< EMaterial< TDomain > > m_spEmMaterial
 Properties of the matherials in the domain. More...
 
SmartPtr< ApproximationSpace< TDomain > > m_spVertApproxSpace
 Approximation space for the potential (vertex-centered) space. More...
 

Detailed Description

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

Projection procedure class.

This class implements the procedure for projection of vector fields in the Nedelec-element representation to the divergence-free space and elimination of the non-zero potentials of conductors.

Member Typedef Documentation

◆ a_vert_cond_type

template<typename TDomain , typename TAlgebra >
typedef AInt ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::a_vert_cond_type
private

Type of the attachment and its accessor the base conductor indices.

◆ aa_vert_cond_type

template<typename TDomain , typename TAlgebra >
typedef Grid::VertexAttachmentAccessor<a_vert_cond_type> ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::aa_vert_cond_type
private

◆ algebra_type

template<typename TDomain , typename TAlgebra >
typedef TAlgebra ug::Electromagnetism::NedelecProject< 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::NedelecProject< TDomain, TAlgebra >::domain_type

Type of Domain.

◆ grid_type

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

Type of Grid:

◆ matrix_type

template<typename TDomain , typename TAlgebra >
typedef TAlgebra::matrix_type ug::Electromagnetism::NedelecProject< 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::NedelecProject< TDomain, TAlgebra >::position_type

position type

◆ pot_gf_type

template<typename TDomain , typename TAlgebra >
typedef GridFunction<TDomain, TPotAlgebra> ug::Electromagnetism::NedelecProject< 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::NedelecProject< 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::NedelecProject< 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::NedelecProject< TDomain, TAlgebra >::subset_handler_type

Type of subset handler.

◆ this_type

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

This type.

◆ TPotAlgebra

template<typename TDomain , typename TAlgebra >
typedef CPUAlgebra ug::Electromagnetism::NedelecProject< 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::NedelecProject< TDomain, TAlgebra >::vector_type

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

Constructor & Destructor Documentation

◆ NedelecProject()

◆ ~NedelecProject()

template<typename TDomain , typename TAlgebra >
ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::~NedelecProject ( )
inline

Member Function Documentation

◆ alloc_DVFs()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::alloc_DVFs ( const TDomain &  domain,
pot_gf_type aux_rhs 
)
private

Allocates memory for the DVFs associated with the ungrounded conductors.

Allocates memory for the DVFs associated with the conductors that do not touch the Dirichlet boundary

Parameters
[in]domainthe domain
[in]aux_rhsgrid function for the auxiliary rhs

References pcl::ProcessCommunicator::allreduce(), ug::GridFunction< class, class >::approx_space(), ug::DoFDistribution::begin(), ug::GridFunction< class, class >::dof_distribution(), ug::DoFDistribution::end(), SmartPtr< class, FreePolicy >::get(), grid(), PCL_RO_LOR, ug::PointerConstArray< class >::size(), and ug::SubsetGroup::size().

◆ apply() [1/3]

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::apply ( GridFunction< TDomain, TAlgebra > &  u,
const char *  fct_names 
)

Performs the projection of given functions.

Performs the projection of given functions

Parameters
[in]uthe grid function to project
[in]fct_namesthe function names

References ug::DoFDistributionInfoProvider::fct_grp_by_name(), and UG_CATCH_THROW.

◆ apply() [2/3]

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::apply ( GridFunction< TDomain, TAlgebra > &  u,
const FunctionGroup fctGrp 
)

◆ apply() [3/3]

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::apply ( vector_type vec)

Performs the projection of all functions in a grid function.

Performs the projection of all functions in a grid function

Parameters
[in]vecthe grid function to project

References ug::FunctionGroup::add_all(), ug::DoFDistributionInfoProvider::dof_distribution_info(), and UG_THROW.

◆ assemble_div()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::assemble_div ( const TDomain &  domain,
const DoFDistribution edgeDD,
const vector_type u,
size_t  fct,
const DoFDistribution vertDD,
pot_vector_type div,
DenseVector< VariableArray1< number > > &  charge 
)
private

Assembles the weak divergence operator on the whole grid.

Assembles the weak divergence operator on the whole grid and computes the charges of the base conductors.

Parameters
[in]domainthe domain
[in]edgeDDthe edge DD
[in]uthe vector to compute div for
[in]fctfunction in u to compute div for
[in]vertDDthe vertex DD
[out]divfor the weak divergence of u
[out]chargecharges of the conductors

References pcl::ProcessCommunicator::allreduce(), PCL_RO_SUM, ug::PST_ADDITIVE, ParallelVector< Vector< double > >::set(), and ParallelVector< Vector< double > >::set_storage_type().

◆ clear_div_in_conductors()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::clear_div_in_conductors ( const DoFDistribution vertDD,
pot_vector_type div,
DenseVector< VariableArray1< number > > &  charge 
)
private

Sets the div operator to 0 in conductors (for all elements of the same type)

Sets the div operator to 0 in conductors for all elements of the same type and computes the charges of the base conductors. The charges of the base conductors are the sums of the divergences in all the interface vertices between the base conductors and the insulators. As the divergence inside of the conductors is set to zero, we sum over all the nodes in the closure of the conductors. Remark: Only full-dimensional elements should be considered.

Parameters
[in]vertDDthe vertex DD
[out]divto update the weak divergence of u
[out]chargecharges of the conductors

References ug::DoFDistribution::algebra_indices(), ug::DoFDistributionInfoProvider::num_subsets(), and UG_THROW.

◆ compute_div()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::compute_div ( SmartPtr< GridFunction< TDomain, TAlgebra > >  sp_u,
const char *  u_fct_name,
SmartPtr< GridFunction< TDomain, TPotAlgebra > >  sp_div 
)

Computes the weak divergence in insulators.

Computes weak divergence in insulators and saves it in a vertex-centered grid function

Parameters
[in]sp_uthe vector field grid function
[in]u_fct_namethe function name of the Nedelec DoFs
[out]sp_divthe grid function for the divergence

References ParallelVector< Vector< double > >::change_storage_type(), ug::GridFunction< class, class >::dof_distribution(), ug::GridFunction< class, class >::domain(), ug::DoFDistributionInfoProvider::fct_id_by_name(), SmartPtr< domain_type >::get(), SmartPtr< class, FreePolicy >::get(), ug::Vector< class >::get(), ug::GridFunction< class, class >::grid_level(), ug::DoFDistributionInfoProvider::local_finite_element_id(), ug::LFEID::NEDELEC, ug::PST_CONSISTENT, ug::LFEID::type(), UG_CATCH_THROW, and UG_THROW.

◆ compute_DVF_potential_coeffs()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::compute_DVF_potential_coeffs ( const TDomain &  domain,
DoFDistribution vertDD 
)
private

Computes the potential coefficients.

Computes the potential coefficients

Parameters
[in]domainthe domain
[in]vertDDthe vertex DD

References pcl::ProcessCommunicator::allreduce(), ug::DoFDistribution::begin(), ug::DoFDistribution::end(), ug::Invert(), ug::DoFDistribution::multi_grid(), PCL_RO_MAX, PCL_RO_SUM, and ug::SetAttachmentValues().

◆ compute_DVFs()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::compute_DVFs ( pot_gf_type aux_rhs)
private

Computes the Dirichlet vector fields (DVFs)

Computes the Dirichlet vector fields

Parameters
aux_rhsgrid function for the auxiliary rhs

References ug::GridFunction< class, class >::dof_distribution(), ug::PST_ADDITIVE, and ug::PST_CONSISTENT.

◆ damp_DVFs()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::damp_DVFs ( pot_vector_type cor,
const DenseVector< VariableArray1< number > > &  charge 
)
private

Damps the Dirichlet vector fields (DVFs)

Damps the Dirichlet vector fields (DVFs)

Parameters
corthe potential correction to update
[in]chargecharges of the conductors

References ug::VecScaleAdd().

◆ distribute_cor()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::distribute_cor ( const TDomain &  domain,
const DoFDistribution edgeDD,
vector_type u,
size_t  fct,
const DoFDistribution vertDD,
const pot_vector_type cor 
)
private

Updates the grid function by the potential correction.

Computes the edge-centered correction from the vertex-centered (potential) one by applying the gradient operator.

Note that the edge-centered correction is considered only in the closure of insulators. There, it is the gradient of the vertex-centered (potential) correction. But in the conductors (as well as on the Dirichlet boundaries), the vertex-centered correction is constant by construction so that its gradient is 0 any way. (Adjacent conductors are concidered as one conductor in the projection, so that the different factors for the charges in different conductors does not violate the constant value of the vertex-centered correction in any separated piece of conductors.) For this, the function computes the gradient over the whole domain.

Parameters
[in]domainthe domain
[in]edgeDDthe edge DD
[out]uthe vector to correct (to update)
[in]fctfunction in u to update
[in]vertDDthe vertex DD
[in]corthe potential correction for u

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

◆ integrate_div_DVF_elem_type()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::integrate_div_DVF_elem_type ( const TDomain &  domain,
const DoFDistribution vertDD,
const aa_vert_cond_type vert_base_cond,
DenseMatrix< VariableArray2< number > > &  C 
)
private

Integrates div E over boundaries of conductors for elements of the same type.

Integration of div E over boundaries of conductors (for one type of elements) to get the capacity matrix. We compute only the lower triangular part of the capacity matrix: This matrix is symmetric.

Parameters
[in]domainthe domain
[in]vertDDthe vertex DD
[in]vert_base_condindices of the base conductors for every vertex
[out]Cthe matrix to update

References ug::DoFDistribution::algebra_indices(), ug::DoFDistributionInfoProvider::num_subsets(), and UG_THROW.

◆ mark_cond_vert_elem_type()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::mark_cond_vert_elem_type ( DoFDistribution vertDD,
aa_vert_cond_type vert_base_cond 
)
private

Sets the base conductor indices for every vertex (for all elements of the same type)

Initializes 'vert_base_cond' with the indices of the base conductors according to the closure of the conductive subsets. Note that the insulators are marked with -2 (and NOT with -1!), whereas the grounded conductors with -1 and the other conductors with their base conductor index. The attachment must be preinitialized with the default value -2 (for the insulators).

Parameters
[in]vertDDthe vertex DD
[out]vert_base_condindices of the base conductors for every vertex

References ug::DoFDistributionInfoProvider::num_subsets().

◆ project_func()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::project_func ( const SmartPtr< TDomain > &  domain,
const SmartPtr< DoFDistribution > &  edgeDD,
vector_type u,
size_t  fct,
const SmartPtr< DoFDistribution > &  vertDD,
pot_gf_type aux_rhs,
pot_gf_type aux_cor 
)
inlineprivate

Projects one function (i.e. performs the main action):

Projects one function, thus performs the main action of the projection. This function assumes that the object is completely initialized, i.e. the solver and the discretization are initialized.

Parameters
[in]domainthe domain
[in]edgeDDthe edge DD
[in]uthe vector where to project
[in]fctfunction in u to compute div for
[in]vertDDthe vertex DD
[in]aux_rhsa grid function for the rhs of the aux. problem
[in]aux_cora grid function for the sol. of the aux. problem

References SmartPtr< class, FreePolicy >::get(), SmartPtr< TDomain >::get(), and ug::PST_CONSISTENT.

◆ set_dampDVFs()

template<typename TDomain , typename TAlgebra >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::set_dampDVFs ( bool  val)
inline

Sets whether to damp Dirichlet vector fields.

References ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::m_bDampDVFs.

◆ set_Dirichlet()

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

Sets the Dirichlet boundary.

Parameters
spDirichletthe Dirichlet BC object

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

◆ weak_div_elem_type()

template<typename TDomain , typename TAlgebra >
template<typename TElem >
void ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::weak_div_elem_type ( const TDomain &  domain,
const DoFDistribution edgeDD,
const vector_type u,
size_t  fct,
const DoFDistribution vertDD,
pot_vector_type div 
)
private

Computes the weak divergence (for all elements of the same type)

Assembles the weak divergence operator for all elements of the same type. Remark: Only full-dimensional elements should be considered.

Parameters
[in]domainthe domain
[in]edgeDDthe edge DD
[in]uthe vector to compute div for
[in]fctfunction in u to compute div for
[in]vertDDthe vertex DD
divto update the weak divergence of u

References ug::DoFDistribution::algebra_indices(), ug::DoFDistribution::dof_indices(), ug::DoFRef(), ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::local_div_matrix(), ug::MatVecMult(), ug::DoFDistributionInfoProvider::num_subsets(), and UG_THROW.

Member Data Documentation

◆ m_auxLaplaceAss

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

Global discretization of the auxiliary equations.

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

◆ m_auxLaplaceOp

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

Matrix of the discretization.

◆ m_auxLaplaceRHS

template<typename TDomain , typename TAlgebra >
SmartPtr<AuxLaplaceRHS> ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::m_auxLaplaceRHS
private

Rhs and BC for the auxiliary problem.

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

◆ m_auxLocLaplace

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

Local discretization of the auxiliary equations.

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

◆ m_bDampDVFs

template<typename TDomain , typename TAlgebra >
bool ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::m_bDampDVFs
private

Whether to damp the Dirichlet vector fields (DVFs)

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

◆ m_DVF_phi

template<typename TDomain , typename TAlgebra >
std::vector<pot_gf_type *> ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::m_DVF_phi
private

◆ m_potCoeff

template<typename TDomain , typename TAlgebra >
DenseMatrix<VariableArray2<number> > ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::m_potCoeff
private

Potential coefficients:

◆ m_potSolver

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

Solver for the auxliliary equations.

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

◆ m_spDirichlet

template<typename TDomain , typename TAlgebra >
SmartPtr<EMDirichlet<TDomain, TAlgebra> > ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::m_spDirichlet
private

◆ m_spEmMaterial

template<typename TDomain , typename TAlgebra >
SmartPtr<EMaterial<TDomain> > ug::Electromagnetism::NedelecProject< TDomain, TAlgebra >::m_spEmMaterial
private

Properties of the matherials in the domain.

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

◆ m_spVertApproxSpace

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

Approximation space for the potential (vertex-centered) space.

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

◆ WDim

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

world dimention


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