Doxygen API reference documentation for ideal.II
|
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... | |
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:
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.
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.
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.
spacetime_vector | The Trilinos vector containing the slab space-time values |
space_vector | The resulting Trilinos vector of spatial values at dof_index |
dof_index | The index of the DoF to be extracted. |
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.
spacetime_vector | The vector containing the slab space-time values |
space_vector | The resulting vector of spatial values at dof_index |
dof_index | The index of the DoF to be extracted. |
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.
dof_handler | The DoFHandler used to calculate the spacetime_vector. |
spacetime_vector | The Trilinos vector containing the slab space-time values. |
space_vector | The resulting Trilinos vector of spatial values at dof_index. |
t | The time point to extract from. |
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.
dof_handler | The DoFHandler used to calculate the spacetime_vector. |
spacetime_vector | The vector containing the slab space-time values. |
space_vector | The resulting vector of spatial values at dof_index. |
t | The time point to extract from. |
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.
space_relevant_dofs | The IndexSet of locally relevant dofs in space. |
dof_handler | The space-time slab dof handler the vector is indexed on |
boundary_component | The boundary id corresponding to the Dirichlet conditions |
boundary_function | The function that describes the Dirichlet data |
spacetime_constraints | The constraints to add the boundary values to |
component_mask | A component mask to only apply boundary values to certain FE components |
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.
dof_handler | The space-time slab dof handler the vector is indexed on |
boundary_component | The boundary id corresponding to the Dirichlet conditions |
boundary_function | The function that describes the Dirichlet data |
spacetime_constraints | The constraints to add the boundary values to |
component_mask | A component mask to only apply boundary values to certain FE components |
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.
dof_handler | The space-time slab dof handler the vector is indexed on |
boundary_component | The boundary id corresponding to the Dirichlet conditions |
boundary_function | The function that describes the Dirichlet data |
spacetime_constraints | The constraints to add the boundary values to |
component_mask | A component mask to only apply boundary values to certain FE components |