Doxygen API reference documentation for ideal.II
Functions
idealii::slab::VectorTools Namespace Reference

Collection of functions working on space-time slab Vectors. More...

Functions

template<int dim, typename Number >
void interpolate_boundary_values (idealii::slab::DoFHandler< dim > &dof_handler, const dealii::types::boundary_id boundary_component, dealii::Function< dim, Number > &boundary_function, std::shared_ptr< dealii::AffineConstraints< Number >> spacetime_constraints, const dealii::ComponentMask &component_mask=dealii::ComponentMask())
 Compute space-time constraints on the solution corresponding to Dirichlet conditions. This function iterates over all temporal degrees of freedom and sets the boundary values of the corresponding subvector to the value of the specified at the time of the temporal dof. More...
 
template<int dim, typename Number = double>
void interpolate_boundary_values (dealii::IndexSet space_relevant_dofs, idealii::slab::DoFHandler< dim > &dof_handler, const dealii::types::boundary_id boundary_component, dealii::Function< dim, Number > &boundary_function, std::shared_ptr< dealii::AffineConstraints< Number >> spacetime_constraints, const dealii::ComponentMask &component_mask=dealii::ComponentMask())
 Compute space-time constraints on the solution corresponding to Dirichlet conditions for the locally relevant set of space dofs. This function iterates over all temporal degrees of freedom and sets the boundary values of the corresponding subvector to the value of the specified at the time of the temporal dof. More...
 
template<int dim, typename Number >
void project_boundary_values_curl_conforming_l2 (idealii::slab::DoFHandler< dim > &dof_handler, unsigned int first_vector_component, dealii::Function< dim, Number > &boundary_function, const dealii::types::boundary_id boundary_component, std::shared_ptr< dealii::AffineConstraints< Number >> spacetime_constraints)
 Compute space-time constraints on the solution corresponding to Hcurl Dirichlet conditions. This function iterates over all temporal degrees of freedom and sets the boundary values of the corresponding subvector to the value of the specified at the time of the temporal dof. More...
 
void extract_subvector_at_time_dof (const dealii::Vector< double > &spacetime_vector, dealii::Vector< double > &space_vector, unsigned int dof_index)
 Get the spatial subvector of a specific temporal dof of the corresponding slab. More...
 
void extract_subvector_at_time_dof (const dealii::TrilinosWrappers::MPI::Vector &spacetime_vector, dealii::TrilinosWrappers::MPI::Vector &space_vector, unsigned int dof_index)
 Get the spatial Trilinos subvector of a specific temporal dof of the corresponding slab. More...
 
template<int dim>
void extract_subvector_at_time_point (slab::DoFHandler< dim > &dof_handler, const dealii::Vector< double > &spacetime_vector, dealii::Vector< double > &space_vector, const double t)
 Get the spatial subvector at a specific time point of the corresponding slab. The result is calculated by linear combination of each temporal dof vector according to the underlying temporal finite element. More...
 
template<int dim>
void extract_subvector_at_time_point (slab::DoFHandler< dim > &dof_handler, const dealii::TrilinosWrappers::MPI::Vector &spacetime_vector, dealii::TrilinosWrappers::MPI::Vector &space_vector, const double t)
 Get the spatial Trilinos subvector at a specific time point of the corresponding slab. The result is calculated by linear combination of each temporal dof vector according to the underlying temporal finite element. More...
 
template<int dim>
double calculate_L2L2_squared_error_on_slab (slab::DoFHandler< dim > &dof_handler, dealii::Vector< double > &spacetime_vector, dealii::Function< dim, double > &exact_solution, spacetime::Quadrature< dim > &quad)
 calculate the L2 inner product of (u-u_{kh}) with itself. More...
 
template<int dim>
double calculate_L2L2_squared_error_on_slab (slab::DoFHandler< dim > &dof_handler, dealii::TrilinosWrappers::MPI::Vector &spacetime_vector, dealii::Function< dim, double > &exact_solution, spacetime::Quadrature< dim > &quad)
 calculate the L2 inner product of (u-u_{kh}) with itself. More...
 

Detailed Description

Collection of functions working on space-time slab Vectors.

These functions provide utilities for manipulating space-time indexed vectors of normally stationary objects provided by deal.II.

Examples are:

Function Documentation

◆ calculate_L2L2_squared_error_on_slab() [1/2]

template<int dim>
double idealii::slab::VectorTools::calculate_L2L2_squared_error_on_slab ( slab::DoFHandler< dim > &  dof_handler,
dealii::TrilinosWrappers::MPI::Vector &  spacetime_vector,
dealii::Function< dim, double > &  exact_solution,
spacetime::Quadrature< dim > &  quad 
)

calculate the L2 inner product of (u-u_{kh}) with itself.

@parameter dof_handler The slab::DoFHandler describing the dof distribution of the space-time vector @parameter spacetime_vector The approximate solution u_{kh} @parameter exact_solution The function describing the exact solution @parameter quad The quadrature formula to use when calculating the integrals.

◆ calculate_L2L2_squared_error_on_slab() [2/2]

template<int dim>
double idealii::slab::VectorTools::calculate_L2L2_squared_error_on_slab ( slab::DoFHandler< dim > &  dof_handler,
dealii::Vector< double > &  spacetime_vector,
dealii::Function< dim, double > &  exact_solution,
spacetime::Quadrature< dim > &  quad 
)

calculate the L2 inner product of (u-u_{kh}) with itself.

@parameter dof_handler The slab::DoFHandler describing the dof distribution of the space-time vector @parameter spacetime_vector The approximate solution u_{kh} @parameter exact_solution The function describing the exact solution @parameter quad The quadrature formula to use when calculating the integrals.

◆ extract_subvector_at_time_dof() [1/2]

void idealii::slab::VectorTools::extract_subvector_at_time_dof ( const dealii::TrilinosWrappers::MPI::Vector &  spacetime_vector,
dealii::TrilinosWrappers::MPI::Vector &  space_vector,
unsigned int  dof_index 
)

Get the spatial Trilinos subvector of a specific temporal dof of the corresponding slab.

Parameters
spacetime_vectorThe Trilinos vector containing the slab space-time values
space_vectorThe resulting Trilinos vector of spatial values at dof_index
dof_indexThe index of the DoF to be extracted.

◆ extract_subvector_at_time_dof() [2/2]

void idealii::slab::VectorTools::extract_subvector_at_time_dof ( const dealii::Vector< double > &  spacetime_vector,
dealii::Vector< double > &  space_vector,
unsigned int  dof_index 
)

Get the spatial subvector of a specific temporal dof of the corresponding slab.

Parameters
spacetime_vectorThe vector containing the slab space-time values
space_vectorThe resulting vector of spatial values at dof_index
dof_indexThe index of the DoF to be extracted.

◆ extract_subvector_at_time_point() [1/2]

template<int dim>
void idealii::slab::VectorTools::extract_subvector_at_time_point ( slab::DoFHandler< dim > &  dof_handler,
const dealii::TrilinosWrappers::MPI::Vector &  spacetime_vector,
dealii::TrilinosWrappers::MPI::Vector &  space_vector,
const double  t 
)

Get the spatial Trilinos subvector at a specific time point of the corresponding slab. The result is calculated by linear combination of each temporal dof vector according to the underlying temporal finite element.

Note
If t is not in the temporal interval of the slab, the resulting vector will be 0.
Parameters
dof_handlerThe DoFHandler used to calculate the spacetime_vector.
spacetime_vectorThe Trilinos vector containing the slab space-time values.
space_vectorThe resulting Trilinos vector of spatial values at dof_index.
tThe time point to extract from.

◆ extract_subvector_at_time_point() [2/2]

template<int dim>
void idealii::slab::VectorTools::extract_subvector_at_time_point ( slab::DoFHandler< dim > &  dof_handler,
const dealii::Vector< double > &  spacetime_vector,
dealii::Vector< double > &  space_vector,
const double  t 
)

Get the spatial subvector at a specific time point of the corresponding slab. The result is calculated by linear combination of each temporal dof vector according to the underlying temporal finite element.

Note
If t is not in the temporal interval of the slab, the resulting vector will be 0.
Parameters
dof_handlerThe DoFHandler used to calculate the spacetime_vector.
spacetime_vectorThe vector containing the slab space-time values.
space_vectorThe resulting vector of spatial values at dof_index.
tThe time point to extract from.

◆ interpolate_boundary_values() [1/2]

template<int dim, typename Number = double>
void idealii::slab::VectorTools::interpolate_boundary_values ( dealii::IndexSet  space_relevant_dofs,
idealii::slab::DoFHandler< dim > &  dof_handler,
const dealii::types::boundary_id  boundary_component,
dealii::Function< dim, Number > &  boundary_function,
std::shared_ptr< dealii::AffineConstraints< Number >>  spacetime_constraints,
const dealii::ComponentMask &  component_mask = dealii::ComponentMask() 
)

Compute space-time constraints on the solution corresponding to Dirichlet conditions for the locally relevant set of space dofs. This function iterates over all temporal degrees of freedom and sets the boundary values of the corresponding subvector to the value of the specified at the time of the temporal dof.

Parameters
space_relevant_dofsThe IndexSet of locally relevant dofs in space.
dof_handlerThe space-time slab dof handler the vector is indexed on
boundary_componentThe boundary id corresponding to the Dirichlet conditions
boundary_functionThe function that describes the Dirichlet data
spacetime_constraintsThe constraints to add the boundary values to
component_maskA component mask to only apply boundary values to certain FE components
Note
The set_time() method of boundary_function is called internally, so if the function is time-dependent use get_time() in your implementation of the function to obtain t.

◆ interpolate_boundary_values() [2/2]

template<int dim, typename Number >
void idealii::slab::VectorTools::interpolate_boundary_values ( idealii::slab::DoFHandler< dim > &  dof_handler,
const dealii::types::boundary_id  boundary_component,
dealii::Function< dim, Number > &  boundary_function,
std::shared_ptr< dealii::AffineConstraints< Number >>  spacetime_constraints,
const dealii::ComponentMask &  component_mask = dealii::ComponentMask() 
)

Compute space-time constraints on the solution corresponding to Dirichlet conditions. This function iterates over all temporal degrees of freedom and sets the boundary values of the corresponding subvector to the value of the specified at the time of the temporal dof.

Parameters
dof_handlerThe space-time slab dof handler the vector is indexed on
boundary_componentThe boundary id corresponding to the Dirichlet conditions
boundary_functionThe function that describes the Dirichlet data
spacetime_constraintsThe constraints to add the boundary values to
component_maskA component mask to only apply boundary values to certain FE components
Note
The set_time() method of boundary_function is called internally, so if the function is time-dependent use get_time() in your implementation of the function to obtain t.

◆ project_boundary_values_curl_conforming_l2()

template<int dim, typename Number >
void idealii::slab::VectorTools::project_boundary_values_curl_conforming_l2 ( idealii::slab::DoFHandler< dim > &  dof_handler,
unsigned int  first_vector_component,
dealii::Function< dim, Number > &  boundary_function,
const dealii::types::boundary_id  boundary_component,
std::shared_ptr< dealii::AffineConstraints< Number >>  spacetime_constraints 
)

Compute space-time constraints on the solution corresponding to Hcurl Dirichlet conditions. This function iterates over all temporal degrees of freedom and sets the boundary values of the corresponding subvector to the value of the specified at the time of the temporal dof.

Parameters
dof_handlerThe space-time slab dof handler the vector is indexed on
boundary_componentThe boundary id corresponding to the Dirichlet conditions
boundary_functionThe function that describes the Dirichlet data
spacetime_constraintsThe constraints to add the boundary values to
component_maskA component mask to only apply boundary values to certain FE components
Note
The set_time() method of boundary_function is called internally, so if the function is time-dependent use get_time() in your implementation of the function to obtain t.