Hermes2D  3.0
Hermes::Hermes2D::RefMap Class Reference

Represents the reference mapping. More...

#include <refmap.h>

+ Inheritance diagram for Hermes::Hermes2D::RefMap:

Public Member Functions

void set_quad_2d (Quad2D *quad_2d)
 
Quad2Dget_quad_2d () const
 Returns the current quadrature points.
 
virtual void set_active_element (Element *e)
 
double3 * get_tangent (int edge, int order=-1)
 
double * get_phys_x (int order)
 
double * get_phys_y (int order)
 
bool is_jacobian_const () const
 
int get_inv_ref_order () const
 Returns the increase in the integration order due to the reference map.
 
double get_const_jacobian () const
 
double2x2 * get_const_inv_ref_map ()
 
double * get_jacobian (int order)
 
double2x2 * get_inv_ref_map (int order)
 
void inv_ref_map_at_point (double xi1, double xi2, double &x, double &y, double2x2 &m)
 Calculates the inverse Jacobi matrix of reference map at a particular point (xi1, xi2).
 
virtual void push_transform (int son)
 See Transformable::push_transform()
 
virtual void pop_transform ()
 See Transformable::pop_transform()
 
void force_transform (uint64_t sub_idx, Trf *ctm)
 For internal use only.
 
- Public Member Functions inherited from Hermes::Hermes2D::Transformable
Elementget_active_element () const
 
virtual void set_transform (uint64_t idx)
 
uint64_t get_transform () const
 

Static Public Member Functions

static void untransform (Element *e, double x, double y, double &xi1, double &xi2)
 
static Elementelement_on_physical_coordinates (bool use_MeshHashGrid, MeshSharedPtr mesh, double x, double y, double *x_reference=nullptr, double *y_reference=nullptr)
 
static bool is_element_on_physical_coordinates (Element *e, double x, double y, double *x_reference=nullptr, double *y_reference=nullptr)
 
static void set_element_iro_cache (Element *element)
 

Additional Inherited Members

- Static Public Attributes inherited from Hermes::Hermes2D::Transformable
static const unsigned int H2D_MAX_TRN_LEVEL = 15
 If this changes, NeighborSearch::H2D_MAX_NEIGHBORS must change too.
 
- Protected Member Functions inherited from Hermes::Hermes2D::Transformable
virtual 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 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 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 inherited from Hermes::Hermes2D::Transformable
static const uint64_t H2D_MAX_IDX = (1ULL << 3 * H2D_MAX_TRN_LEVEL) - 1
 

Detailed Description

Represents the reference mapping.

RefMap represents the mapping from the reference to the physical element. Its main task is to provide both the (inverse) reference mapping matrix and its jacobian, precalculated in quadrature points. Other functions include the calculation of integration points positions in the physical domain and the calculation of edge tangents in 1D integration points.

Definition at line 40 of file refmap.h.

Member Function Documentation

Element * Hermes::Hermes2D::RefMap::element_on_physical_coordinates ( bool  use_MeshHashGrid,
MeshSharedPtr  mesh,
double  x,
double  y,
double *  x_reference = nullptr,
double *  y_reference = nullptr 
)
static

Returns the element pointer located at physical coordinates x, y.

Parameters
[in]xPhysical x-coordinate.
[in]yPhysical y-coordinate.
[in]x_referenceOptional parameter, in which the x-coordinate of x in the reference domain will be returned.
[in]y_referenceOptional parameter, in which the y-coordinate of y in the reference domain will be returned.

Definition at line 765 of file refmap.cpp.

double2x2* Hermes::Hermes2D::RefMap::get_const_inv_ref_map ( )
inline

If the reference map is constant, this is the fast way to obtain its inverse matrix.

Definition at line 115 of file refmap.h.

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

double Hermes::Hermes2D::RefMap::get_const_jacobian ( ) const
inline

If the jacobian of the reference map is constant, this is the fast way to obtain it.

Definition at line 108 of file refmap.h.

double2x2* Hermes::Hermes2D::RefMap::get_inv_ref_map ( int  order)
inline

Returns the inverse matrices of the reference map precalculated at the integration points of the specified order. Intended for non-constant jacobian elements.

Definition at line 134 of file refmap.h.

Referenced by Hermes::Hermes2D::init_fn_preallocated(), and Hermes::Hermes2D::Mesh::initial_single_check().

double* Hermes::Hermes2D::RefMap::get_jacobian ( int  order)
inline

Returns the jacobian of the reference map precalculated at the integration points of the specified order. Intended for non-constant jacobian elements.

Definition at line 122 of file refmap.h.

Referenced by Hermes::Hermes2D::Mesh::initial_single_check().

double * Hermes::Hermes2D::RefMap::get_phys_x ( int  order)

Returns the x-coordinates of the integration points transformed to the physical domain of the element. Intended for integrals containing spatial variables.

Definition at line 38 of file refmap.cpp.

Referenced by Hermes::Hermes2D::init_geom_surf_allocated(), Hermes::Hermes2D::init_geom_vol_allocated(), Hermes::Hermes2D::VonMisesFilter::precalculate(), and Hermes::Hermes2D::Views::Orderizer::process_space().

double * Hermes::Hermes2D::RefMap::get_phys_y ( int  order)

Returns he y-coordinates of the integration points transformed to the physical domain of the element. Intended for integrals containing spatial variables.

Definition at line 45 of file refmap.cpp.

Referenced by Hermes::Hermes2D::init_geom_surf_allocated(), Hermes::Hermes2D::init_geom_vol_allocated(), and Hermes::Hermes2D::Views::Orderizer::process_space().

double3 * Hermes::Hermes2D::RefMap::get_tangent ( int  edge,
int  order = -1 
)

Returns the triples[x, y, norm] of the tangent to the specified (possibly curved) edge at the 1D integration points along the edge. The maximum 1D quadrature rule is used by default, but the user may specify his own order. In this case, the edge pseudo-order is expected (as returned by Quad2D::get_edge_points).

Definition at line 52 of file refmap.cpp.

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

bool Hermes::Hermes2D::RefMap::is_element_on_physical_coordinates ( Element e,
double  x,
double  y,
double *  x_reference = nullptr,
double *  y_reference = nullptr 
)
static

Find out if the coordinatex [x,y] lie in the element e.

Parameters
[in]xPhysical x-coordinate.
[in]yPhysical y-coordinate.
[in]x_referenceOptional parameter, in which the x-coordinate of x in the reference domain will be returned.
[in]y_referenceOptional parameter, in which the y-coordinate of y in the reference domain will be returned.

Definition at line 670 of file refmap.cpp.

Referenced by Hermes::Hermes2D::MeshHashGridElement::getElement().

bool Hermes::Hermes2D::RefMap::is_jacobian_const ( ) const
inline

Returns true if the jacobian of the reference map is constant (which is the case for non-curvilinear triangular elements), false otherwise.

Definition at line 95 of file refmap.h.

Referenced by Hermes::Hermes2D::init_fn_preallocated(), and Hermes::Hermes2D::Mesh::initial_single_check().

void Hermes::Hermes2D::RefMap::set_active_element ( Element e)
virtual
void Hermes::Hermes2D::RefMap::set_quad_2d ( Quad2D quad_2d)

Sets the quadrature points in which the reference map will be evaluated.

Parameters
quad_2d[in]The quadrature points.

Definition at line 63 of file refmap.cpp.

Referenced by Hermes::Hermes2D::PostProcessing::VolumetricIntegralCalculator< Scalar >::calculate(), Hermes::Hermes2D::PostProcessing::SurfaceIntegralCalculator< Scalar >::calculate(), and Hermes::Hermes2D::Views::Orderizer::process_space().

void Hermes::Hermes2D::RefMap::untransform ( Element e,
double  x,
double  y,
double &  xi1,
double &  xi2 
)
static

Transforms physical coordinates x, y from the element e back to the reference domain. If the point (x, y) does not lie in e, then (xi1, xi2) will not lie in the reference domain.

Definition at line 542 of file refmap.cpp.


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