Plugins
ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem > Class Template Reference

Tool kit for the Whitney-1 (Nedelec) based FE discretization of the rot-rot operators. More...

#include <nedelec_local_ass.h>

Public Types

typedef TDomain::grid_type grid_type
 type of the grid More...
 
typedef TDomain::position_type position_type
 type of the geometric positions (WDim-vectors) More...
 
typedef reference_element_traits< TElem >::reference_element_type ref_elem_type
 type of reference element More...
 
typedef LagrangeP1< ref_elem_typeW0_shapes_type
 shapes and derivatives of the Whitney-0 (Lagrange P1) shape functions More...
 

Static Public Member Functions

static void curl (const TDomain *domain, TElem *elem, const position_type *corners, const number dofs[], MathVector< WDim > &curl_vec)
 computes the curl of the grid functions (in 2d represented as \((z, 0)\) instead of \((0, 0, z)\)) More...
 
static void get_shapes (const TDomain *domain, TElem *elem, const position_type *corners, const MathVector< dim > local, MathVector< WDim > shapes[])
 computes the Nedelec shapes at a given point More...
 
static void interpolate (const TDomain *domain, TElem *elem, const position_type *corners, const number dofs[], const MathVector< dim > local[], const size_t n_pnt, MathVector< WDim > values[])
 computes of the values of the grid functions More...
 
static void local_div_matrix (const TDomain *domain, TElem *elem, const position_type *corners, number B[][numEdges])
 assembles the discrete weak div operator More...
 
static void local_mass (const TDomain *domain, TElem *elem, const position_type *corners, number M[maxNumEdges][maxNumEdges])
 assembles the local mass matrix of the Nedelec elements More...
 
static void local_stiffness_and_mass (const TDomain *domain, TElem *elem, const position_type *corners, number S[maxNumEdges][maxNumEdges], number M[maxNumEdges][maxNumEdges])
 assembles the local stiffness and the mass matrices of the rot-rot operator More...
 

Static Public Attributes

static const int dim = ref_elem_type::dim
 dimension of reference element More...
 
static const size_t maxNumEdges = (size_t) element_list_traits<typename domain_traits<WDim>::DimElemList>::maxEdges
 max. number of the edges of the full-dimensional elements in the domain More...
 
static const size_t numCorners = ref_elem_type::numCorners
 total number of the corners More...
 
static const size_t numEdges = ref_elem_type::numEdges
 total number of the edges More...
 
static const int WDim = TDomain::dim
 world dimention More...
 

Detailed Description

template<typename TDomain, typename TElem>
class ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >

Tool kit for the Whitney-1 (Nedelec) based FE discretization of the rot-rot operators.

Class for the local discretization of the rot-rot operator using the Nedelec-type-1 (Whitney-1) elements. This class may not contain any member except for static functions, so no instances of this class need to be created.

This basic class template does not contain any numerical algorithms as they are specific for different types of the elements. Cf. the specializations for simplices (

See also
NedelecT1_LDisc_forSimplex).
Remarks
Note that the definition of the Whitney-1 function depends on the ordering of the nodes in the element. More precisely, for a different ordering, the shape function \(w^{(1)}_e\) may have a different sign. Thus, the ordering should be consistent in all the elements: If two elements \(T_1\) and \(T_2\) share the edge \(e\) with the ends at corners \(i\) and \(j\) (that are shared, too) then this edge should be \(e = (i, j)\) for both the elements (and not for example \(e = (i, j)\) for \(T_1\) but \(e = (j, i)\) for \(T_2\)). In all the other respects, the ordering does not play any essential role. In particular, the singes of the degrees of freedom in the solution depend on the global ordering of the corners (but the absolute values of these dofs are invariant). For this, the solution of the discretized system may be only considered in connection with the given ordering of the corners (in the edges).

Member Typedef Documentation

◆ grid_type

template<typename TDomain , typename TElem >
typedef TDomain::grid_type ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::grid_type

type of the grid

◆ position_type

template<typename TDomain , typename TElem >
typedef TDomain::position_type ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::position_type

type of the geometric positions (WDim-vectors)

◆ ref_elem_type

template<typename TDomain , typename TElem >
typedef reference_element_traits<TElem>::reference_element_type ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::ref_elem_type

type of reference element

◆ W0_shapes_type

template<typename TDomain , typename TElem >
typedef LagrangeP1<ref_elem_type> ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::W0_shapes_type

shapes and derivatives of the Whitney-0 (Lagrange P1) shape functions

Member Function Documentation

◆ curl()

template<typename TDomain , typename TElem >
static void ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::curl ( const TDomain *  domain,
TElem *  elem,
const position_type corners,
const number  dofs[],
MathVector< WDim > &  curl_vec 
)
inlinestatic

computes the curl of the grid functions (in 2d represented as \((z, 0)\) instead of \((0, 0, z)\))

Parameters
[in]domainthe domain
[in]elemelement
[in]cornersarray of the global corner coordinates
[in]dofsarrays of values of the Nedelec degrees of freedom
[out]curl_vecwhere to store the computed curl

References UG_THROW.

Referenced by ug::Electromagnetism::NedelecInterpolation< TDomain, 2, 2 >::curl(), ug::Electromagnetism::NedelecInterpolation< TDomain, 3, 3 >::curl(), and ug::Electromagnetism::NedelecInterpolation< TDomain, 3, 3 >::curl_flux().

◆ get_shapes()

template<typename TDomain , typename TElem >
static void ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::get_shapes ( const TDomain *  domain,
TElem *  elem,
const position_type corners,
const MathVector< dim local,
MathVector< WDim shapes[] 
)
inlinestatic

computes the Nedelec shapes at a given point

Parameters
[in]domainthe domain
[in]elemelement
[in]cornersarray of the global corner coordinates
[in]locallocal coordinates of the point where to compute
[out]shapesarray for the shapes

References UG_THROW.

Referenced by ug::Electromagnetism::NedelecProlongationMatrixHelper< TDomain, TAlgebra, TElem >::assemble_prolongation_matrix().

◆ interpolate()

template<typename TDomain , typename TElem >
static void ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::interpolate ( const TDomain *  domain,
TElem *  elem,
const position_type corners,
const number  dofs[],
const MathVector< dim local[],
const size_t  n_pnt,
MathVector< WDim values[] 
)
inlinestatic

computes of the values of the grid functions

Parameters
[in]domainthe domain
[in]elemelement
[in]cornersarray of the global corner coordinates
[in]dofsarrays of values of the Nedelec degrees of freedom
[in]locallocal coordinates of the points where to compute
[in]n_pntnumber of the points where to compute
[out]valueswhere to store the computed n_pnt values

References UG_THROW.

Referenced by ug::Electromagnetism::NedelecInterpolation< TDomain, 2, 2 >::value(), and ug::Electromagnetism::NedelecInterpolation< TDomain, 3, 3 >::value().

◆ local_div_matrix()

template<typename TDomain , typename TElem >
static void ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::local_div_matrix ( const TDomain *  domain,
TElem *  elem,
const position_type corners,
number  B[][numEdges] 
)
inlinestatic

assembles the discrete weak div operator

Parameters
[in]domainthe domain
[in]elemelement
[in]cornersarray of the global corner coordinates
[out]Blocal weak divergence operator matrix

References UG_THROW.

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

◆ local_mass()

template<typename TDomain , typename TElem >
static void ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::local_mass ( const TDomain *  domain,
TElem *  elem,
const position_type corners,
number  M[maxNumEdges][maxNumEdges] 
)
inlinestatic

assembles the local mass matrix of the Nedelec elements

Parameters
[in]domainthe domain
[in]elemelement
[in]cornersarray of the global corner coordinates
[out]Mlocal mass matrix

References UG_THROW.

Referenced by ug::Electromagnetism::CalcVolPowerElemHelperClass< TGridFunc, TElem >::calc_elem_power().

◆ local_stiffness_and_mass()

template<typename TDomain , typename TElem >
static void ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::local_stiffness_and_mass ( const TDomain *  domain,
TElem *  elem,
const position_type corners,
number  S[maxNumEdges][maxNumEdges],
number  M[maxNumEdges][maxNumEdges] 
)
inlinestatic

assembles the local stiffness and the mass matrices of the rot-rot operator

Parameters
[in]domainthe domain
[in]elemelement
[in]cornersarray of the global corner coordinates
[out]Slocal stiffness matrix
[out]Mlocal mass matrix

References UG_THROW.

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

Member Data Documentation

◆ dim

template<typename TDomain , typename TElem >
const int ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::dim = ref_elem_type::dim
static

dimension of reference element

◆ maxNumEdges

template<typename TDomain , typename TElem >
const size_t ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::maxNumEdges = (size_t) element_list_traits<typename domain_traits<WDim>::DimElemList>::maxEdges
static

max. number of the edges of the full-dimensional elements in the domain

◆ numCorners

template<typename TDomain , typename TElem >
const size_t ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::numCorners = ref_elem_type::numCorners
static

total number of the corners

◆ numEdges

template<typename TDomain , typename TElem >
const size_t ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::numEdges = ref_elem_type::numEdges
static

total number of the edges

◆ WDim

template<typename TDomain , typename TElem >
const int ug::Electromagnetism::NedelecT1_LDisc< TDomain, TElem >::WDim = TDomain::dim
static

world dimention


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