Doxygen API reference documentation for ideal.II
Public Member Functions | Public Attributes | List of all members
idealii::spacetime::FEFaceValues< dim > Class Template Reference

Evaluation of the tensor-product space-time basis functions on spatial element faces. More...

#include <spacetime_fe_values.hh>

Public Member Functions

 FEFaceValues (DG_FiniteElement< dim > &fe, Quadrature< dim - 1 > &quad, const dealii::UpdateFlags uflags, const dealii::UpdateFlags additional_flags)
 Constructor of the FEValues class. More...
 
void reinit_space (const typename dealii::TriaIterator< dealii::DoFCellAccessor< dim, dim, false >> &cell_space, const unsigned int face_no)
 Reinitialize all objects of the underlying spatial FEValues object. This function calls reinit(cell_space) of the spatial FEValues object. More...
 
void reinit_time (const typename dealii::TriaIterator< dealii::DoFCellAccessor< 1, 1, false >> &cell_time)
 Reinitialize all objects of the underlying temporal FEValues object. This function calls reinit(cell_time) of the temporal FEValues object. More...
 
double shape_value (unsigned int function_no, unsigned int point_no)
 Value of the space-time shape function at spacetime-quadrature point. More...
 
dealii::FEValuesViews::Scalar< dim >::value_type scalar_value (const typename dealii::FEValuesExtractors::Scalar &extractor, unsigned int function_no, unsigned int point_no)
 Value of the space-time shape function of a scalar finite element component. More...
 
dealii::FEValuesViews::Vector< dim >::value_type vector_value (const typename dealii::FEValuesExtractors::Vector &extractor, unsigned int function_no, unsigned int point_no)
 Value of the space-time shape function of a vector-valued finite element component. More...
 
double time_quadrature_point (unsigned int quadrature_point)
 Get the temporal quadrature point of the given space-time quadrature index. More...
 
dealii::Point< dim > space_quadrature_point (unsigned int quadrature_point)
 Get the spatial quadrature point of the given space-time quadrature index. More...
 
const dealii::Tensor< 1, dim > & space_normal_vector (unsigned int i)
 Get the normal vector at the spatial face. More...
 
double JxW (const unsigned int quadrature_point)
 Mapped space-time quadrature weight. More...
 
void get_local_dof_indices (std::vector< dealii::types::global_dof_index > &indices)
 Local space-time DoF indices of the current space-time element. More...
 
std::shared_ptr< dealii::FEFaceValues< dim > > spatial ()
 The underlying spatial FEValues object. More...
 
std::shared_ptr< dealii::FEValues< 1 > > temporal ()
 The underlying temporal FEValues object. More...
 

Public Attributes

unsigned int n_quadrature_points
 Number of space-time quadrature points per element.
 

Detailed Description

template<int dim>
class idealii::spacetime::FEFaceValues< dim >

Evaluation of the tensor-product space-time basis functions on spatial element faces.

This class supplies common derivatives of the space-time basis functions by multiplying the corresponding spatial and temporal basis functions in the given quadrature points.

In practice spatial derivatives are handled by the underlying spatial dealii::FEFaceValues object and temporal derivatives are handled by the underlying temporal dealii::FEValues.

Constructor & Destructor Documentation

◆ FEFaceValues()

template<int dim>
idealii::spacetime::FEFaceValues< dim >::FEFaceValues ( DG_FiniteElement< dim > &  fe,
Quadrature< dim - 1 > &  quad,
const dealii::UpdateFlags  uflags,
const dealii::UpdateFlags  additional_flags 
)

Constructor of the FEValues class.

Parameters
feThe underlying space-time finite element description class.
quadThe space-time quadrature formula to be used.
uflagsThe update flags to be used during the reinit calls.
additional_flagsAdditional update flags that are only appended to the FEFace Object (e.g.: update_normal_vectors)

Member Function Documentation

◆ get_local_dof_indices()

template<int dim>
void idealii::spacetime::FEFaceValues< dim >::get_local_dof_indices ( std::vector< dealii::types::global_dof_index > &  indices)

Local space-time DoF indices of the current space-time element.

Parameters
Avector of indices to save the result to.

◆ JxW()

template<int dim>
double idealii::spacetime::FEFaceValues< dim >::JxW ( const unsigned int  quadrature_point)

Mapped space-time quadrature weight.

Parameters
quadrature_pointspace-time quadrature index.

◆ reinit_space()

template<int dim>
void idealii::spacetime::FEFaceValues< dim >::reinit_space ( const typename dealii::TriaIterator< dealii::DoFCellAccessor< dim, dim, false >> &  cell_space,
const unsigned int  face_no 
)

Reinitialize all objects of the underlying spatial FEValues object. This function calls reinit(cell_space) of the spatial FEValues object.

Parameters
cell_spaceIterator pointing to the current element in space.

◆ reinit_time()

template<int dim>
void idealii::spacetime::FEFaceValues< dim >::reinit_time ( const typename dealii::TriaIterator< dealii::DoFCellAccessor< 1, 1, false >> &  cell_time)

Reinitialize all objects of the underlying temporal FEValues object. This function calls reinit(cell_time) of the temporal FEValues object.

Parameters
cell_timeIterator pointing to the current element in time.

◆ scalar_value()

template<int dim>
dealii::FEValuesViews::Scalar<dim>::value_type idealii::spacetime::FEFaceValues< dim >::scalar_value ( const typename dealii::FEValuesExtractors::Scalar &  extractor,
unsigned int  function_no,
unsigned int  point_no 
)

Value of the space-time shape function of a scalar finite element component.

This function passes the extractor to the underlying spatial FEValues object and then calls the value function of the resulting view.

Parameters
extractorScalar extractor defining the scalar finite element component to be evaluated.
function_noThe number of the space-time function/dof to be evaluated.
point_noThe number of the quadrature point to evaluate at.

◆ shape_value()

template<int dim>
double idealii::spacetime::FEFaceValues< dim >::shape_value ( unsigned int  function_no,
unsigned int  point_no 
)

Value of the space-time shape function at spacetime-quadrature point.

Parameters
function_noThe number of the space-time function/dof to be evaluated.
point_noThe number of the quadrature point to evaluate at.

◆ space_normal_vector()

template<int dim>
const dealii::Tensor<1, dim>& idealii::spacetime::FEFaceValues< dim >::space_normal_vector ( unsigned int  i)

Get the normal vector at the spatial face.

Parameters
iThe space-time index of the quadrature point

◆ space_quadrature_point()

template<int dim>
dealii::Point<dim> idealii::spacetime::FEFaceValues< dim >::space_quadrature_point ( unsigned int  quadrature_point)

Get the spatial quadrature point of the given space-time quadrature index.

Parameters
quadrature_pointspace-time quadrature index.

◆ spatial()

template<int dim>
std::shared_ptr<dealii::FEFaceValues<dim> > idealii::spacetime::FEFaceValues< dim >::spatial ( )

The underlying spatial FEValues object.

Returns
A shared pointer to the spatial FEValues object.

◆ temporal()

template<int dim>
std::shared_ptr<dealii::FEValues<1> > idealii::spacetime::FEFaceValues< dim >::temporal ( )

The underlying temporal FEValues object.

Returns
A shared pointer to the temporal FEValues object.

◆ time_quadrature_point()

template<int dim>
double idealii::spacetime::FEFaceValues< dim >::time_quadrature_point ( unsigned int  quadrature_point)

Get the temporal quadrature point of the given space-time quadrature index.

Parameters
quadrature_pointspace-time quadrature index.

◆ vector_value()

template<int dim>
dealii::FEValuesViews::Vector<dim>::value_type idealii::spacetime::FEFaceValues< dim >::vector_value ( const typename dealii::FEValuesExtractors::Vector &  extractor,
unsigned int  function_no,
unsigned int  point_no 
)

Value of the space-time shape function of a vector-valued finite element component.

This function passes the extractor to the underlying spatial FEValues object and then calls the value function of the resulting view.

Parameters
extractorScalar extractor defining the scalar finite element component to be evaluated.
function_noThe number of the space-time function/dof to be evaluated.
point_noThe number of the quadrature point to evaluate at.

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