Hermes2D
2.0
|
Represents the reference mapping. More...
#include <refmap.h>
Public Member Functions | |
void | set_quad_2d (Quad2D *quad_2d) |
Quad2D * | get_quad_2d () const |
Returns the current quadrature points. | |
const Quad1D * | get_quad_1d () const |
Returns the 1D quadrature for use in surface integrals. | |
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 |
double | get_const_jacobian () const |
double * | get_jacobian (int order) |
int | get_inv_ref_order () const |
Returns the increase in the integration order due to the reference map. | |
double2x2 * | get_inv_ref_map (int order) |
double2x2 * | get_const_inv_ref_map () |
![]() | |
Element * | get_active_element () const |
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 Element * | element_on_physical_coordinates (const Mesh *mesh, double x, double y, double *x_reference=NULL, double *y_reference=NULL) |
Public Attributes | |
H1ShapesetJacobi | ref_map_shapeset |
PrecalcShapeset | ref_map_pss |
Friends | |
template<typename T > | |
class | MeshFunction |
template<typename T > | |
class | DiscreteProblem |
template<typename T > | |
class | DiscreteProblemLinear |
template<typename T > | |
class | Solution |
template<typename T > | |
class | ExactSolution |
template<typename T > | |
class | ExactSolutionScalar |
template<typename T > | |
class | ExactSolutionVector |
template<typename T > | |
class | Adapt |
template<typename T > | |
class | KellyTypeAdapt |
class | Views::Orderizer |
class | Views::Vectorizer |
class | Views::Linearizer |
template<typename T > | |
class | Global |
class | VonMisesFilter |
template<typename T > | |
class | Func |
template<typename T > | |
class | Geom |
Geom< double > * | init_geom_vol (RefMap *rm, const int order) |
Init element geometry for volumetric integrals. | |
Geom< double > * | init_geom_surf (RefMap *rm, SurfPos *surf_pos, const int order) |
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). | |
Additional Inherited Members | |
![]() | |
void | reset_transform () |
Empties the stack, loads identity transform. | |
double | get_transform_jacobian () const |
Trf * | get_ctm () const |
unsigned int | get_depth () const |
![]() | |
static void | push_transforms (std::set< Transformable * > &transformables, int son) |
static void | pop_transforms (std::set< Transformable * > &transformables) |
![]() | |
Element * | element |
The active element. | |
Trf * | ctm |
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 const unsigned int | H2D_MAX_TRN_LEVEL = 15 |
static const uint64_t | H2D_MAX_IDX = (1ULL << 3 * H2D_MAX_TRN_LEVEL) - 1 |
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.
|
static |
Returns the element pointer located at physical coordinates x, y.
[in] | x | Physical x-coordinate. |
[in] | y | Physical y-coordinate. |
[in] | x_reference | Optional parameter, in which the x-coordinate of x in the reference domain will be returned. |
[in] | y_reference | Optional parameter, in which the y-coordinate of y in the reference domain will be returned. |
Definition at line 711 of file refmap.cpp.
double2x2 * Hermes::Hermes2D::RefMap::get_const_inv_ref_map | ( | ) |
If the reference map is constant, this is the fast way to obtain its inverse matrix.
Definition at line 69 of file refmap.cpp.
Referenced by Hermes::Hermes2D::init_fn().
double Hermes::Hermes2D::RefMap::get_const_jacobian | ( | ) | const |
If the jacobian of the reference map is constant, this is the fast way to obtain it.
Definition at line 62 of file refmap.cpp.
Referenced by Hermes::Hermes2D::DiscreteProblem< Scalar >::init_geometry_points().
double2x2 * Hermes::Hermes2D::RefMap::get_inv_ref_map | ( | int | order | ) |
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 88 of file refmap.cpp.
Referenced by Hermes::Hermes2D::init_fn().
double * Hermes::Hermes2D::RefMap::get_jacobian | ( | int | order | ) |
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 76 of file refmap.cpp.
Referenced by Hermes::Hermes2D::Adapt< Scalar >::eval_error(), Hermes::Hermes2D::Adapt< Scalar >::eval_error_norm(), Hermes::Hermes2D::KellyTypeAdapt< Scalar >::eval_volumetric_estimator(), and Hermes::Hermes2D::DiscreteProblem< Scalar >::init_geometry_points().
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 109 of file refmap.cpp.
Referenced by Hermes::Hermes2D::init_geom_surf(), Hermes::Hermes2D::init_geom_vol(), and Hermes::Hermes2D::VonMisesFilter::precalculate().
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 120 of file refmap.cpp.
Referenced by Hermes::Hermes2D::init_geom_surf(), and Hermes::Hermes2D::init_geom_vol().
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 133 of file refmap.cpp.
Referenced by Hermes::Hermes2D::init_geom_surf().
bool Hermes::Hermes2D::RefMap::is_jacobian_const | ( | ) | const |
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 49 of file refmap.cpp.
Referenced by Hermes::Hermes2D::init_fn(), and Hermes::Hermes2D::DiscreteProblem< Scalar >::init_geometry_points().
|
virtual |
Initializes the reference map for the specified element. Must be called prior to using all other functions in the class.
Reimplemented from Hermes::Hermes2D::Transformable.
Definition at line 158 of file refmap.cpp.
Referenced by Hermes::Hermes2D::DiscreteProblem< Scalar >::assemble_one_state().
void Hermes::Hermes2D::RefMap::set_quad_2d | ( | Quad2D * | quad_2d | ) |
Sets the quadrature points in which the reference map will be evaluated.
quad_2d[in] | The quadrature points. |
Definition at line 151 of file refmap.cpp.
Referenced by Hermes::Hermes2D::DiscreteProblem< Scalar >::assemble_one_DG_state().
|
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 585 of file refmap.cpp.
Referenced by element_on_physical_coordinates().