Plugins
ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim > Class Template Reference

Interpolation of the Nedelec dofs and their curls. More...

#include <nedelec_local_ass.h>

Public Types

typedef TDomain::position_type position_type
 type of the geometric positions (WDim-vectors) More...
 

Static Public Member Functions

static void curl (const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], MathVector< WDim > &curl_vec)
 computes curl of the function More...
 
static number curl_flux (const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< WDim > &normal, const position_type pnt, number &flux)
 computes flux of the curl through a given plane in an element More...
 
static void value (const TDomain *domain, GridObject *elem, const position_type corners[], const number dofs[], const MathVector< refDim > local[], const size_t n_pnt, MathVector< WDim > values[])
 computes the values at given points More...
 

Static Public Attributes

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...
 

Detailed Description

template<typename TDomain, int refDim, int WDim = TDomain::dim>
class ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim >

Interpolation of the Nedelec dofs and their curls.

This class implements the transformation of the Nedelec element into the vector field.

Template Parameters
TDomaintype of the domain
refDimdimensionality of the reference element
WDimdimensionality of the domain (Do not specify it yourself! It is only for the specializations.)

Member Typedef Documentation

◆ position_type

template<typename TDomain , int refDim, int WDim = TDomain::dim>
typedef TDomain::position_type ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim >::position_type

type of the geometric positions (WDim-vectors)

Member Function Documentation

◆ curl()

template<typename TDomain , int refDim, int WDim = TDomain::dim>
static void ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim >::curl ( const TDomain *  domain,
GridObject elem,
const position_type  corners[],
const number  dofs[],
MathVector< WDim > &  curl_vec 
)
inlinestatic

computes curl of the function

This function computes the value of the curl operator for the Nedelec representation. Curl is constant over elements. Note that the result is represented as a vector, not as a Whitney-2-form. In 2d, where the result has always the form (0, 0, z), it is represented as (z, 0).

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 n_pnt values

References ug::GridObject::reference_object_id(), and UG_THROW.

Referenced by ug::Electromagnetism::EddyCurrentCurlEDependentCmpUserData< TImpl, ReIm, TGFunc >::eval_and_deriv(), and ug::Electromagnetism::NedelecCurlData< TGridFunc >::eval_and_deriv().

◆ curl_flux()

template<typename TDomain , int refDim, int WDim = TDomain::dim>
static number ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim >::curl_flux ( const TDomain *  domain,
GridObject elem,
const position_type  corners[],
const number  dofs[],
const MathVector< WDim > &  normal,
const position_type  pnt,
number flux 
)
inlinestatic

computes flux of the curl through a given plane in an element

This function computes the flux of the curl through a given plane inside of a given grid element. The plane is identified by a point on it and the normal. The flux is multiplied by the norm of the given normal vector (i.e. specify the unit normal to get the standard flux). The function returns the area of the intersection if the element is intersected by the plane. Otherwise the function returns exactly 0.0.

Parameters
[in]domainthe domain
[in]elemelement
[in]cornersarray of the global corner coordinates
[in]dofsarrays of values of the Nedelec degrees of freedom
[in]normalnormal to the plane
[in]pntpoint on the plane (identifying the plane)
[out]fluxthe flux

References ug::GridObject::reference_object_id(), and UG_THROW.

Referenced by ug::Electromagnetism::calc_EMF().

◆ value()

template<typename TDomain , int refDim, int WDim = TDomain::dim>
static void ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim >::value ( const TDomain *  domain,
GridObject elem,
const position_type  corners[],
const number  dofs[],
const MathVector< refDim >  local[],
const size_t  n_pnt,
MathVector< WDim >  values[] 
)
inlinestatic

computes the values at given points

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::GridObject::reference_object_id(), and UG_THROW.

Referenced by ug::Electromagnetism::ComputeElemFluxHelper< TGridFunc, TElem >::compute_elem_flux(), ug::Electromagnetism::NedelecGridFunctionData< TGridFunc >::eval_and_deriv(), ug::Electromagnetism::NedelecSigmaEData< TGridFunc >::eval_and_deriv(), and ug::Electromagnetism::EddyCurrentHeat< TGridFunc >::eval_and_deriv().

Member Data Documentation

◆ maxNumEdges

template<typename TDomain , int refDim, int WDim = TDomain::dim>
const size_t ug::Electromagnetism::NedelecInterpolation< TDomain, refDim, WDim >::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


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