Doxygen API reference documentation for ideal.II
|
Collection of functions working on degrees of freedom. More...
Functions | |
template<int dim> | |
void | make_upwind_sparsity_pattern (DoFHandler< dim > &dof, dealii::DynamicSparsityPattern &st_dsp, std::shared_ptr< dealii::AffineConstraints< double >> space_constraints=std::make_shared< dealii::AffineConstraints< double >>(), const bool keep_constrained_dofs=true) |
Construction of a sparsity pattern with lower off diagonal jump terms. More... | |
template<int dim> | |
void | make_upwind_sparsity_pattern (DoFHandler< dim > &dof, const dealii::Table< 2, dealii::DoFTools::Coupling > &space_couplings, dealii::DynamicSparsityPattern &st_dsp, std::shared_ptr< dealii::AffineConstraints< double >> space_constraints=std::make_shared< dealii::AffineConstraints< double >>(), const bool keep_constrained_dofs=true) |
Construction of a sparsity pattern with lower off diagonal jump terms. More... | |
template<int dim> | |
void | make_downwind_sparsity_pattern (DoFHandler< dim > &dof, dealii::DynamicSparsityPattern &st_dsp, std::shared_ptr< dealii::AffineConstraints< double >> space_constraints=std::make_shared< dealii::AffineConstraints< double >>(), const bool keep_constrained_dofs=true) |
Construction of a sparsity pattern with upper off diagonal jump terms. More... | |
template<int dim> | |
void | make_downwind_sparsity_pattern (DoFHandler< dim > &dof, const dealii::Table< 2, dealii::DoFTools::Coupling > &space_couplings, dealii::DynamicSparsityPattern &st_dsp, std::shared_ptr< dealii::AffineConstraints< double >> space_constraints=std::make_shared< dealii::AffineConstraints< double >>(), const bool keep_constrained_dofs=true) |
Construction of a sparsity pattern with upper off diagonal jump terms. More... | |
template<int dim, typename Number > | |
void | make_hanging_node_constraints (idealii::slab::DoFHandler< dim > &dof_handler, std::shared_ptr< dealii::AffineConstraints< Number >> spacetime_constraints) |
Construction of space-time hanging node constraints. More... | |
template<int dim, typename Number > | |
void | make_hanging_node_constraints (dealii::IndexSet &space_relevant_dofs, idealii::slab::DoFHandler< dim > &dof_handler, std::shared_ptr< dealii::AffineConstraints< Number >> spacetime_constraints) |
Construction of parallel distributed space-time hanging node constraints. More... | |
template<int dim> | |
dealii::IndexSet | extract_locally_relevant_dofs (slab::DoFHandler< dim > &dof_handler) |
Extract global space-time dof indices active on the current DoFHandler. More... | |
Collection of functions working on degrees of freedom.
These functions provide utilities for calculating spacetime variants of normally stationary objects provided by deal.II.
In detail this means spreading stationary information over all temporal degrees of freedom in this slab. Examples are:
dealii::IndexSet idealii::slab::DoFTools::extract_locally_relevant_dofs | ( | slab::DoFHandler< dim > & | dof_handler | ) |
Extract global space-time dof indices active on the current DoFHandler.
For DoFHandlers with a parallel triangulation this function returns the locally owned dofs and all dofs belonging to other processors that are located at a locally owned cell.s
void idealii::slab::DoFTools::make_downwind_sparsity_pattern | ( | DoFHandler< dim > & | dof, |
const dealii::Table< 2, dealii::DoFTools::Coupling > & | space_couplings, | ||
dealii::DynamicSparsityPattern & | st_dsp, | ||
std::shared_ptr< dealii::AffineConstraints< double >> | space_constraints = std::make_shared<dealii::AffineConstraints<double>>() , |
||
const bool | keep_constrained_dofs = true |
||
) |
Construction of a sparsity pattern with upper off diagonal jump terms.
This functions works like the previous one but the spatial pattern is constructed using the supplied couplings.
space_couplings | In coupled problems some components don't couple in the weak formulation and consequently don't needs entries in the sparsity pattern. |
An example would be the Stokes system where pressure ansatz function is not multiplied by the pressure test function, which leads to an empty diagonal block.
void idealii::slab::DoFTools::make_downwind_sparsity_pattern | ( | DoFHandler< dim > & | dof, |
dealii::DynamicSparsityPattern & | st_dsp, | ||
std::shared_ptr< dealii::AffineConstraints< double >> | space_constraints = std::make_shared<dealii::AffineConstraints<double>>() , |
||
const bool | keep_constrained_dofs = true |
||
) |
Construction of a sparsity pattern with upper off diagonal jump terms.
This functions constructs the tensor product between a spatial sparsity pattern and a temporal upwind sparsity.
The spatial pattern is constructed by calling the function make_sparsity_pattern() with the spatial DoFHandler of dof
, the space_contraints
and keep_constrained_dofs
.
The temporal pattern is initially constructed as a block-diagonal sparsity pattern by the corresponding deal.II method for the temporal DoFHandler of dof
. If the number of temporal elements in this slab is greater than one, additional off diagonal entries for the jump terms between a temporal element and its direct right/later neighbor are added.
For Gauss-Lobatto support points in time this is only a coupling between the final temporal dof of the current element and the initial temporal dof of the neighbor.
For Gauss-Legendre support points in time this is the complete first upper off diagonal block.
dof | A shared pointer to the DoFHandler object to base the patterns on. |
st_dsp | The dynamic sparsity pattern that will afterwards include the spacetime pattern. |
space_constraints | The spatial constraints handed to the spatial pattern function |
keep_constrained_dofs | When using distribute_local_to_global() of the spacetime constraints off diagonal entries of constrained dofs will not be written. So it is possible to not include these in the sparsity pattern. For more details see the functions in the dealii::DoFTools namespace. |
void idealii::slab::DoFTools::make_hanging_node_constraints | ( | dealii::IndexSet & | space_relevant_dofs, |
idealii::slab::DoFHandler< dim > & | dof_handler, | ||
std::shared_ptr< dealii::AffineConstraints< Number >> | spacetime_constraints | ||
) |
Construction of parallel distributed space-time hanging node constraints.
This function internally constructs spatial hanging node constraints based on the spatial part of @dof_handler.
The spacetime constraints are then obtained by offsetting the constraints by the total number of dofs in the spatial dof handler.
dof_handler | |
spacetime_constraints |
void idealii::slab::DoFTools::make_hanging_node_constraints | ( | idealii::slab::DoFHandler< dim > & | dof_handler, |
std::shared_ptr< dealii::AffineConstraints< Number >> | spacetime_constraints | ||
) |
Construction of space-time hanging node constraints.
This function internally constructs spatial hanging node constraints based on the spatial part of @dof_handler.
The spacetime constraints are then obtained by offsetting the constraints by the total number of dofs in the spatial dof handler.
dof_handler | |
spacetime_constraints |
void idealii::slab::DoFTools::make_upwind_sparsity_pattern | ( | DoFHandler< dim > & | dof, |
const dealii::Table< 2, dealii::DoFTools::Coupling > & | space_couplings, | ||
dealii::DynamicSparsityPattern & | st_dsp, | ||
std::shared_ptr< dealii::AffineConstraints< double >> | space_constraints = std::make_shared<dealii::AffineConstraints<double>>() , |
||
const bool | keep_constrained_dofs = true |
||
) |
Construction of a sparsity pattern with lower off diagonal jump terms.
This functions works like the previous one butthe spatial pattern is constructed using the supplied couplings.
space_couplings | In coupled problems some components don't couple in the weak formulation and consequently don't needs entries in the sparsity pattern. |
An example would be the Stokes system where pressure ansatz function is not multiplied by the pressure test function, which leads to an empty diagonal block.
void idealii::slab::DoFTools::make_upwind_sparsity_pattern | ( | DoFHandler< dim > & | dof, |
dealii::DynamicSparsityPattern & | st_dsp, | ||
std::shared_ptr< dealii::AffineConstraints< double >> | space_constraints = std::make_shared<dealii::AffineConstraints<double>>() , |
||
const bool | keep_constrained_dofs = true |
||
) |
Construction of a sparsity pattern with lower off diagonal jump terms.
This functions constructs the tensor product between a spatial sparsity pattern and a temporal upwind sparsity.
The spatial pattern is constructed by calling the function make_sparsity_pattern() with the spatial DoFHandler of dof
, the space_contraints
and keep_constrained_dofs
.
The temporal pattern is initially constructed as a block-diagonal sparsity pattern by the corresponding deal.II method for the temporal DoFHandler of dof
. If the number of temporal elements in this slab is greater than one, additional off diagonal entries for the jump terms between a temporal element and its direct left/earlier neighbor are added.
For Gauss-Lobatto support points in time this is only a coupling between the initial temporal dof of the current element and the final temporal dof of the neighbor.
For Gauss-Legendre support points in time this is the complete first lower off diagonal block.
dof | A shared pointer to the DoFHandler object to base the patterns on. |
st_dsp | The dynamic sparsity pattern that will afterwards include the spacetime pattern. |
space_constraints | The spatial constraints handed to the spatial pattern function |
keep_constrained_dofs | When using distribute_local_to_global() of the spacetime constraints off diagonal entries of constrained dofs will not be written. So it is possible to not include these in the sparsity pattern. For more details see the functions in the dealii::DoFTools namespace. |