hexed 0.3.0
 
Loading...
Searching...
No Matches
hexed::Element Class Reference

Stores data associated with one mesh element. More...

#include <Element.hpp>

Inheritance diagram for hexed::Element:
hexed::Kernel_element hexed::Deformed_element

Public Member Functions

 Element (Storage_params, std::vector< int > pos={}, double mesh_size=1., int ref_level=0, Mat<> origin_arg=Mat<>::Zero(3), int aniso_ref_level=0)
 
virtual bool get_is_deformed ()
 
 Element (const Element &)=delete
 Can't copy an Element. Doing so would have to either duplicate or break vertex connections, both of which seem error prone.
 
Elementoperator= (const Element &)=delete
 
Storage_params storage_params ()
 
virtual std::vector< double > position (const Basis &, int i_qpoint)
 
std::vector< double > face_position (const Basis &, int i_face, int i_face_qpoint)
 obtains face position based on interior qpoint positions (as defined by position())
 
virtual void set_jacobian (const Basis &basis)
 
double nominal_size () const override
 nominal edge length of the element before any vertex motion
 
int refinement_level ()
 indicates how many times this element has been isotropically refined
 
int aniso_ref_level ()
 indicates how many times this element has been anisotropically refined
 
std::vector< int > nominal_position ()
 
double * stage (int i_stage)
 pointer to state data for i_stageth Runge-Kutta stage.
 
double * advection_state ()
 
double * time_step_scale () override
 pointer to scaling factor for local time step.
 
double * bulk_av_coef ()
 layout: [i_qpoint]
 
double * laplacian_av_coef ()
 layout: [i_qpoint]
 
double * art_visc_forcing ()
 layout: [i_forcing][i_qpoint]
 
virtual double * node_adjustments ()
 overriden by Deformed_element
 
int mask () const override
 returns whether the element is included in the masked mesh.
 
virtual double jacobian (int i_dim, int j_dim, int i_qpoint)
 Compute the Jacobian matrix.
 
virtual double jacobian_determinant (int i_qpoint)
 determinant of jacobian
 
Vertexvertex (int i_vertex)
 
template<int i_dim>
double & vertex_position (int i_vertex)
 
void push_shareable_value (std::function< double(Element &, int i_vertex)>)
 functions to communicate with nodal neighbors
 
void fetch_shareable_value (std::function< double &(Element &, int i_vertex)> access_fun, std::function< double(Mat<>)> reduction=Vertex::vector_max)
 
double & vertex_time_step_scale (int i_vertex) override
 Time step scale at the vertices. TSS in the interior is set by interpolating this.
 
double & vertex_elwise_av (int i_vertex)
 
double & vertex_fix_admis_coef (int i_vertex)
 
void set_needs_smooth (bool)
 sets the Vertex::Transferable_ptr::needs_smooth of the vertices
 
void set_face (int i_face, double *data)
 
bool is_connected (int i_face)
 
double * state () override
 pointer to the data where the state variables are stored
 
double * residual_cache () override
 
double * face (int i_face, bool is_ldg) override
 
bool deformed () const override
 where the extrapolated face data is stored
 
double * reference_level_normals () override
 whether this element is deformed
 
double * jacobian_determinant () override
 Jacobian determinant of reference -> physical coordinate transform.
 
double * kernel_face_normal (int i_face) override
 where the extrapolated face area-weighted normal vectors are stored
 
double & uncert () override
 place to store some uncertainty metric
 

Public Attributes

std::array< int, 6 > face_record
 for algorithms to book-keep information related to faces
 
double uncertainty = 0
 Pointer to state data at faces. Must be populated by user.
 
Mutual_ptr< Element, Treetree
 Tree this element was created from
 
bool unrefinement_locked = false
 if this is set to true, Mesh_interface::update() won't unrefine it
 
bool snapping_problem = false
 if true, this element has a face on the surface which was not properly snapped
 
bool needs_snapping = true
 once any faces of this element have been snapped to the surface, set this to false
 
const Mat origin
 origin which integer coordinates are relative to
 
Lock lock
 for any tasks where multiple threads might access an element simultaneously
 
- Public Attributes inherited from hexed::Kernel_element
int record = 0
 for algorithms to book-keep general information
 

Static Public Attributes

static constexpr bool is_deformed = false
 is this Element subclass deformed?
 

Protected Member Functions

 Element (Storage_params, std::vector< int > pos, double mesh_size, int ref_level, Mat<> origin_arg, bool mobile_vertices, int aniso_r_level)
 

Protected Attributes

Storage_params params
 
int n_dim
 
std::vector< int > _nom_pos
 
double _nom_sz
 
int _r_level
 
int _aniso_r_level
 
std::vector< Vertex::Transferable_ptrvertices
 

Detailed Description

Stores data associated with one mesh element.

Container only—does not have implementations of or information about the basis and algorithms. This class represents a Cartesian (i.e., regular) element. See also derived class Deformed_element.

Constructor & Destructor Documentation

◆ Element()

hexed::Element::Element ( Storage_params params_arg,
std::vector< int > pos = {},
double mesh_size = 1.,
int ref_level = 0,
Mat<> origin_arg = Mat<>::Zero(3),
int aniso_ref_level = 0 )

The Storage_params defines the amount of storage that must be allocated. pos specifies the position of vertex 0 relative to origin_arg in intervals of the nominal size. The nominal size is defined to be mesh_size/(2^ref_level). The vertices will be spaced at intervals of the nominal size. Only the first n_dim elements of origin_arg are considered.

Member Function Documentation

◆ advection_state()

double * hexed::Element::advection_state ( )

layout: [i_node][i_qpoint]

Note
0 <= i_node < row_size

◆ deformed()

bool hexed::Element::deformed ( ) const
overridevirtual

where the extrapolated face data is stored

Implements hexed::Kernel_element.

◆ face()

double * hexed::Element::face ( int i_face,
bool is_ldg )
overridevirtual

Implements hexed::Kernel_element.

◆ face_position()

std::vector< double > hexed::Element::face_position ( const Basis & basis,
int i_face,
int i_face_qpoint )

obtains face position based on interior qpoint positions (as defined by position())

Note
ignores vertex positions

◆ get_is_deformed()

virtual bool hexed::Element::get_is_deformed ( )
inlinevirtual

for determining whether a pointer is deformed

Reimplemented in hexed::Deformed_element.

◆ jacobian()

double hexed::Element::jacobian ( int i_dim,
int j_dim,
int i_qpoint )
virtual

Compute the Jacobian matrix.

I.e., derivative of i_dimth physical coordinate wrt j_dimth reference coordinate. Trivial for this class, may be non-trivial for derived (see Deformed_element). For convenience, not performance. For high-performance, use double* Deformed_element::jacobian()

Reimplemented in hexed::Deformed_element.

◆ jacobian_determinant() [1/2]

double * hexed::Element::jacobian_determinant ( )
overridevirtual

Jacobian determinant of reference -> physical coordinate transform.

for Cartesian elements, returns nullptr

Implements hexed::Kernel_element.

◆ jacobian_determinant() [2/2]

virtual double hexed::Element::jacobian_determinant ( int i_qpoint)
inlinevirtual

determinant of jacobian

Reimplemented in hexed::Deformed_element.

◆ kernel_face_normal()

double * hexed::Element::kernel_face_normal ( int i_face)
overridevirtual

where the extrapolated face area-weighted normal vectors are stored

for Cartesian elements, returns nullptr

Implements hexed::Kernel_element.

◆ mask()

int hexed::Element::mask ( ) const
inlineoverridevirtual

returns whether the element is included in the masked mesh.

value can be set with Accessible_mesh::set_mask

Implements hexed::Kernel_element.

◆ node_adjustments()

virtual double * hexed::Element::node_adjustments ( )
inlinevirtual

overriden by Deformed_element

Reimplemented in hexed::Deformed_element.

◆ nominal_size()

double hexed::Element::nominal_size ( ) const
inlineoverridevirtual

nominal edge length of the element before any vertex motion

Implements hexed::Kernel_element.

◆ position()

std::vector< double > hexed::Element::position ( const Basis & basis,
int i_qpoint )
virtual

Reimplemented in hexed::Deformed_element.

◆ push_shareable_value()

void hexed::Element::push_shareable_value ( std::function< double(Element &, int i_vertex)> fun)

functions to communicate with nodal neighbors

writes shareable value to vertices so that shared value can be determined

◆ reference_level_normals()

double * hexed::Element::reference_level_normals ( )
overridevirtual

whether this element is deformed

the j_dimth component (in physical space) of the normal vector of the level surface of the i_dimth reference coordinate which passes through the i_qpointth quadrature point. the magnitude of the normal vector is weighted by the surface area in physical space. equivalent definition (note i_dim, j_dim transposed):

\[ J^{-1}_{j_{dim}, i_{dim}} |J| \]

where is jacobian of transformation from reference to physical coordinates at i_qpointth quadrature point. For Cartesian elements, returns nullptr.

layout: [i_dim][j_dim][i_qpoint]

Implements hexed::Kernel_element.

◆ residual_cache()

double * hexed::Element::residual_cache ( )
overridevirtual

Implements hexed::Kernel_element.

◆ set_jacobian()

void hexed::Element::set_jacobian ( const Basis & basis)
virtual

Reimplemented in hexed::Deformed_element.

◆ stage()

double * hexed::Element::stage ( int i_stage)

pointer to state data for i_stageth Runge-Kutta stage.

layout: [i_var][i_qpoint]

◆ state()

double * hexed::Element::state ( )
overridevirtual

pointer to the data where the state variables are stored

includes any non-conservation variables such as artificial viscosity coefficient

Implements hexed::Kernel_element.

◆ time_step_scale()

double * hexed::Element::time_step_scale ( )
overridevirtual

pointer to scaling factor for local time step.

layout: [i_qpoint]

Implements hexed::Kernel_element.

◆ uncert()

double & hexed::Element::uncert ( )
overridevirtual

place to store some uncertainty metric

Implements hexed::Kernel_element.

◆ vertex_time_step_scale()

double & hexed::Element::vertex_time_step_scale ( int i_vertex)
overridevirtual

Time step scale at the vertices. TSS in the interior is set by interpolating this.

Implements hexed::Kernel_element.

Member Data Documentation

◆ uncertainty

double hexed::Element::uncertainty = 0

Pointer to state data at faces. Must be populated by user.

refinement algorithms should set this value to some uncertainty metric


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