Hermes2D  2.0
Hermes::Hermes2D::Function< Scalar > Class Template Referenceabstract

Represents an arbitrary function defined on an element. More...

#include <function.h>

+ Inheritance diagram for Hermes::Hermes2D::Function< Scalar >:

Classes

struct  Node
 

Public Member Functions

 Function ()
 Default constructor.
 
virtual ~Function ()
 
int get_num_components () const
 Returns the number of components of the function being represented by the class.
 
Scalar * get_fn_values (int component=0)
 Returns function values. More...
 
Scalar * get_dx_values (int component=0)
 Returns the x partial derivative. More...
 
Scalar * get_dy_values (int component=0)
 Returns the y partial derivative. More...
 
void get_dx_dy_values (Scalar *&dx, Scalar *&dy, int component=0)
 Returns both x and y partial derivatives. This function provides the both often-used dx and dy values in one call. More...
 
Scalar * get_dxx_values (int component=0)
 Returns the second x partial derivative. More...
 
Scalar * get_dyy_values (int component=0)
 Returns the second y partial derivative. More...
 
Scalar * get_dxy_values (int component=0)
 Returns the second mixed derivative. More...
 
Quad2Dget_quad_2d () const
 Returns the current quadrature points.
 
void set_quad_order (unsigned int order, int mask=H2D_FN_DEFAULT)
 
Scalar * get_values (int a, int b)
 
int get_fn_order () const
 Returns the polynomial degree of the function being represented by the class.
 
- Public Member Functions inherited from Hermes::Hermes2D::Transformable
Elementget_active_element () const
 
void set_transform (uint64_t idx)
 
uint64_t get_transform () const
 
virtual void push_transform (int son)
 

Protected Member Functions

virtual void set_quad_2d (Quad2D *quad_2d)
 Selects the quadrature points in which the function will be evaluated. More...
 
virtual void free ()=0
 Frees all precalculated tables.
 
virtual int get_edge_fn_order (int edge) const
 Returns the polynomial degree of the function at given edge. To be overridden in derived classes. More...
 
virtual void precalculate (int order, int mask)=0
 precalculates the current function at the current integration points.
 
void update_nodes_ptr ()
 
void force_transform (uint64_t sub_idx, Trf *ctm)
 For internal use only.
 
Nodenew_node (int mask, int num_points)
 allocates a new Node structure
 
virtual void handle_overflow_idx ()=0
 
void replace_cur_node (Node *node)
 
- Protected Member Functions inherited from Hermes::Hermes2D::Transformable
virtual void set_active_element (Element *e)
 
virtual void pop_transform ()
 
void reset_transform ()
 Empties the stack, loads identity transform.
 
double get_transform_jacobian () const
 
Trfget_ctm () const
 
unsigned int get_depth () const
 

Static Protected Member Functions

static void check_params (int component, Node *cur_node, int num_components)
 
static void check_table (int component, Node *cur_node, int n, const char *msg)
 
- Static Protected Member Functions inherited from Hermes::Hermes2D::Transformable
static void push_transforms (std::set< Transformable * > &transformables, int son)
 
static void pop_transforms (std::set< Transformable * > &transformables)
 

Protected Attributes

int order
 current function polynomial order
 
int num_components
 number of vector components
 
std::map< uint64_t, LightArray
< Node * > * > * 
sub_tables
 Table of Node tables, for each possible transformation there can be a different Node table.
 
LightArray< Node * > * nodes
 Table of nodes.
 
Nodecur_node
 Current Node.
 
LightArray< Node * > * overflow_nodes
 Nodes for the overflow sub-element transformation.
 
Quad2Dquads [8]
 list of available quadratures
 
int cur_quad
 active quadrature (index into 'quads')
 
int total_mem
 total memory in bytes used by the tables
 
int max_mem
 peak memory usage
 
- Protected Attributes inherited from Hermes::Hermes2D::Transformable
Elementelement
 The active element.
 
Trfctm
 Current sub-element transformation matrix.
 
uint64_t sub_idx
 Sub-element transformation index.
 
Trf stack [21]
 Transformation matrix stack.
 
unsigned int top
 Stack top.
 

Static Protected Attributes

static int idx2mask [6][2]
 index to mask table More...
 
- Static Protected Attributes inherited from Hermes::Hermes2D::Transformable
static const unsigned int H2D_MAX_TRN_LEVEL = 15
 
static const uint64_t H2D_MAX_IDX = (1ULL << 3 * H2D_MAX_TRN_LEVEL) - 1
 

Friends

template<typename T >
class KellyTypeAdapt
 
template<typename T >
class RefinementSelectors::H1ProjBasedSelector
 
template<typename T >
class RefinementSelectors::L2ProjBasedSelector
 
template<typename T >
class RefinementSelectors::HcurlProjBasedSelector
 
template<typename T >
class Adapt
 
class Views::Orderizer
 
class Views::Vectorizer
 
class Views::Linearizer
 
template<typename T >
class DiscontinuousFunc
 
template<typename T >
class DiscreteProblem
 
template<typename T >
class DiscreteProblemLinear
 
template<typename T >
class Global
 
class CurvMap
 
template<typename T >
class Func
 
template<typename T >
class Geom
 
template<typename T >
class Filter
 
template<typename T >
class SimpleFilter
 
template<typename T >
class DXDYFilter
 
class ComplexFilter
 
class VonMisesFilter
 
HERMES_API Geom< double > * init_geom_vol (RefMap *rm, const int order)
 Init element geometry for volumetric integrals.
 
HERMES_API Geom< double > * init_geom_surf (RefMap *rm, SurfPos *surf_pos, const int order)
 
HERMES_API Func< double > * init_fn (PrecalcShapeset *fu, RefMap *rm, const int order)
 Init the shape function for the evaluation of the volumetric/surface integral (transformation of values).
 
template<typename T >
HERMES_API Func< T > * init_fn (MeshFunction< T > *fu, const int order)
 

Detailed Description

template<typename Scalar>
class Hermes::Hermes2D::Function< Scalar >

Represents an arbitrary function defined on an element.

The Function class is an abstraction of a function defined in integration points on an element. You first specify what quadrature tables you want to use (set_quad_2d()) and select an element (Transformable::set_active_element()). Then you select concrete integration points (set_quad_order()) and obtain the function values by calling one of the functions get_fn_values(), get_dx_values(), etc.

This class is a template for RealFunction and ScalarFunction, depending of which type the function values are. For example, shape functions are always Real (see PrecalcShapeset), while the solution can be complex (see Solution).

The design goal for this class is to define a single common interface for functions used as integrands in the weak formulation. It should not matter whether you are integrating a shape function or, for example, a previous solution of the PDE in time-dependent problems. Ideally, you should also be able to apply the bilinear form not only to shape functions during assembling, but also to the solution when calculating energy norms etc. The last feature is unfortunately limited to Real code, because a PDE solution can be complex (hence Solution inherits from ScalarFunction), but shape functions are Real and for efficiency the bilinear form only takes RealFunction arguments.

Since this class inherits from Transformable, you can obtain function values in integration points transformed to sub-areas of the current element (see push_transform(), pop_transform()).

Definition at line 98 of file function.h.

Constructor & Destructor Documentation

template<typename Scalar >
Hermes::Hermes2D::Function< Scalar >::~Function ( )
virtual

Default destructor. Added to delete Function::sub_tables.

Definition at line 53 of file function.cpp.

Member Function Documentation

template<typename Scalar>
void Hermes::Hermes2D::Function< Scalar >::get_dx_dy_values ( Scalar *&  dx,
Scalar *&  dy,
int  component = 0 
)

Returns both x and y partial derivatives. This function provides the both often-used dx and dy values in one call.

Parameters
dx[out]Variable which receives the pointer to the first partial derivatives by x
dy[out]Variable which receives the pointer to the first partial derivatives by y
component[in]The component of the function (0 or 1).

Definition at line 116 of file function.cpp.

Referenced by Hermes::Hermes2D::VonMisesFilter::precalculate().

template<typename Scalar >
Scalar * Hermes::Hermes2D::Function< Scalar >::get_dx_values ( int  component = 0)

Returns the x partial derivative.

Parameters
component[in]The component of the function (0 or 1).
Returns
The x partial derivative of the function at all points of the current integration rule.

Definition at line 102 of file function.cpp.

Referenced by Hermes::Hermes2D::init_fn(), Hermes::Hermes2D::RefinementSelectors::H1ProjBasedSelector< Scalar >::precalc_ref_solution(), and Hermes::Hermes2D::RefinementSelectors::HcurlProjBasedSelector< Scalar >::precalc_ref_solution().

template<typename Scalar >
Scalar * Hermes::Hermes2D::Function< Scalar >::get_dxx_values ( int  component = 0)

Returns the second x partial derivative.

Parameters
component[in]The component of the function (0 or 1).
Returns
The x second partial derivative of the function at all points of the current integration rule.

Definition at line 124 of file function.cpp.

Referenced by Hermes::Hermes2D::init_fn().

template<typename Scalar >
Scalar * Hermes::Hermes2D::Function< Scalar >::get_dxy_values ( int  component = 0)

Returns the second mixed derivative.

Parameters
component[in]The component of the function (0 or 1).
Returns
The second mixed derivative of the function at all points of the current integration rule.

Definition at line 138 of file function.cpp.

Referenced by Hermes::Hermes2D::init_fn().

template<typename Scalar >
Scalar * Hermes::Hermes2D::Function< Scalar >::get_dy_values ( int  component = 0)

Returns the y partial derivative.

Parameters
component[in]The component of the function (0 or 1).
Returns
The y partial derivative of the function at all points of the current integration rule.

Definition at line 109 of file function.cpp.

Referenced by Hermes::Hermes2D::init_fn(), Hermes::Hermes2D::RefinementSelectors::H1ProjBasedSelector< Scalar >::precalc_ref_solution(), and Hermes::Hermes2D::RefinementSelectors::HcurlProjBasedSelector< Scalar >::precalc_ref_solution().

template<typename Scalar >
Scalar * Hermes::Hermes2D::Function< Scalar >::get_dyy_values ( int  component = 0)

Returns the second y partial derivative.

Parameters
component[in]The component of the function (0 or 1).
Returns
The y second partial derivative of the function at all points of the current integration rule.

Definition at line 131 of file function.cpp.

Referenced by Hermes::Hermes2D::init_fn().

template<typename Scalar >
int Hermes::Hermes2D::Function< Scalar >::get_edge_fn_order ( int  edge) const
protectedvirtual

Returns the polynomial degree of the function at given edge. To be overridden in derived classes.

Parameters
edge[in]Edge at which the order should be evaluated. (0-3)

Definition at line 64 of file function.cpp.

template<typename Scalar >
Scalar * Hermes::Hermes2D::Function< Scalar >::get_fn_values ( int  component = 0)
template<typename Scalar >
void Hermes::Hermes2D::Function< Scalar >::set_quad_2d ( Quad2D quad_2d)
protectedvirtual

Selects the quadrature points in which the function will be evaluated.

It is possible to switch back and forth between different quadrature points: no precalculated values are freed. The standard quadrature is always selected by default already.

Parameters
quad_2d[in]The quadrature points.

Reimplemented in Hermes::Hermes2D::ComplexFilter, Hermes::Hermes2D::MeshFunction< Scalar >, Hermes::Hermes2D::MeshFunction< double >, Hermes::Hermes2D::MeshFunction< std::complex< double > >, Hermes::Hermes2D::Filter< Scalar >, Hermes::Hermes2D::Filter< double >, and Hermes::Hermes2D::Filter< std::complex< double > >.

Definition at line 151 of file function.cpp.

Referenced by Hermes::Hermes2D::MeshFunction< Scalar >::set_quad_2d().

template<typename Scalar >
void Hermes::Hermes2D::Function< Scalar >::set_quad_order ( unsigned int  order,
int  mask = H2D_FN_DEFAULT 
)

Activates an integration rule of the specified order. Subsequent calls to get_values(), get_dx_values() etc. will be returning function values at these points.

Parameters
order[in]Integration rule order.
mask[in]A combination of one or more of the constants H2D_FN_VAL, H2D_FN_DX, H2D_FN_DY, H2D_FN_DXX, H2D_FN_DYY, H2D_FN_DXY specifying the values which should be precalculated. The default is H2D_FN_VAL | H2D_FN_DX | H2D_FN_DY. You can also use H2D_FN_ALL to precalculate everything.

Definition at line 76 of file function.cpp.

Referenced by Hermes::Hermes2D::RefinementSelectors::ProjBasedSelector< Scalar >::calc_projection_errors(), Hermes::Hermes2D::init_fn(), and Hermes::Hermes2D::VonMisesFilter::precalculate().

template<typename Scalar >
void Hermes::Hermes2D::Function< Scalar >::update_nodes_ptr ( )
protected

With changed sub-element mapping, there comes the need for a change of the current Node table nodes.

Definition at line 217 of file function.cpp.

Referenced by Hermes::Hermes2D::MeshFunction< Scalar >::pop_transform(), and Hermes::Hermes2D::MeshFunction< Scalar >::push_transform().

Member Data Documentation

template<typename Scalar>
int Hermes::Hermes2D::Function< Scalar >::idx2mask
staticprotected
Initial value:
=
{
{ H2D_FN_VAL_0, H2D_FN_VAL_1 }, { H2D_FN_DX_0, H2D_FN_DX_1 }, { H2D_FN_DY_0, H2D_FN_DY_1 },
{ H2D_FN_DXX_0, H2D_FN_DXX_1 }, { H2D_FN_DYY_0, H2D_FN_DYY_1 }, { H2D_FN_DXY_0, H2D_FN_DXY_1 }
}

index to mask table

Definition at line 240 of file function.h.


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