hexed 0.3.0
 
Loading...
Searching...
No Matches
hexed::Kernel_element Class Referenceabstract

Represents an element as the kernels see it. More...

#include <Kernel_element.hpp>

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

Public Member Functions

virtual int mask () const =0
 
virtual double * state ()=0
 pointer to the data where the state variables are stored
 
virtual double * residual_cache ()=0
 
virtual double * time_step_scale ()=0
 where the convective residual is stored for 2-stage time integration
 
virtual double & vertex_time_step_scale (int i_vertex)=0
 
virtual double nominal_size () const =0
 nominal edge length of the element before any vertex motion
 
virtual double * face (int i_face, bool is_ldg)=0
 
virtual bool deformed () const =0
 where the extrapolated face data is stored
 
virtual double * reference_level_normals ()=0
 whether this element is deformed
 
virtual double * jacobian_determinant ()=0
 Jacobian determinant of reference -> physical coordinate transform.
 
virtual double * kernel_face_normal (int i_face)=0
 where the extrapolated face area-weighted normal vectors are stored
 
virtual double & uncert ()=0
 place to store some uncertainty metric
 

Public Attributes

int record = 0
 for algorithms to book-keep general information
 

Detailed Description

Represents an element as the kernels see it.

Includes only the bare minimum of information needed to run the spatial discretization kernels. This minimizes the number of times where an irrelevant change necessitates recompiling the kernels, which is expensive and thus very annoying. It also helps to better define what the kernels do and what exactly their inputs and outputs are. The kernel is responsible for knowing what the Storage_params and storage order of the element are.

Member Function Documentation

◆ deformed()

virtual bool hexed::Kernel_element::deformed ( ) const
pure virtual

where the extrapolated face data is stored

Implemented in hexed::Deformed_element, hexed::Element, and hexed::Element_new.

◆ jacobian_determinant()

virtual double * hexed::Kernel_element::jacobian_determinant ( )
pure virtual

Jacobian determinant of reference -> physical coordinate transform.

for Cartesian elements, returns nullptr

Implemented in hexed::Deformed_element, hexed::Element, and hexed::Element_new.

◆ kernel_face_normal()

virtual double * hexed::Kernel_element::kernel_face_normal ( int i_face)
pure virtual

where the extrapolated face area-weighted normal vectors are stored

for Cartesian elements, returns nullptr

Implemented in hexed::Deformed_element, hexed::Element, and hexed::Element_new.

◆ mask()

virtual int hexed::Kernel_element::mask ( ) const
pure virtual

Implemented in hexed::Element.

◆ nominal_size()

virtual double hexed::Kernel_element::nominal_size ( ) const
pure virtual

nominal edge length of the element before any vertex motion

Implemented in hexed::Element, and hexed::Element_new.

◆ reference_level_normals()

virtual double * hexed::Kernel_element::reference_level_normals ( )
pure virtual

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]

Implemented in hexed::Deformed_element, hexed::Element, and hexed::Element_new.

◆ state()

virtual double * hexed::Kernel_element::state ( )
pure virtual

pointer to the data where the state variables are stored

includes any non-conservation variables such as artificial viscosity coefficient

Implemented in hexed::Element, and hexed::Element_new.

◆ time_step_scale()

virtual double * hexed::Kernel_element::time_step_scale ( )
pure virtual

where the convective residual is stored for 2-stage time integration

storage for local time step

layout: [i_qpoint]

Implemented in hexed::Element, and hexed::Element_new.

◆ uncert()

virtual double & hexed::Kernel_element::uncert ( )
pure virtual

place to store some uncertainty metric

Implemented in hexed::Element, and hexed::Element_new.

◆ vertex_time_step_scale()

virtual double & hexed::Kernel_element::vertex_time_step_scale ( int i_vertex)
pure virtual
Todo
the kernel should not need this

Implemented in hexed::Element, and hexed::Element_new.


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