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

Actual DoFHandler for a specific slab. More...

#include <slab_dof_handler.hh>

Public Member Functions

 DoFHandler (Triangulation< dim > &tria)
 Constructor linking a slab::Triangulation. More...
 
 DoFHandler (slab::parallel::distributed::Triangulation< dim > &tria)
 Constructor linking a parallel::distributed::slab::Triangulation. More...
 
 DoFHandler (const DoFHandler< dim > &other)
 (shallow) copy constructor. Only the index set and fe support type are actually copied. The underlying pointers will point to the same dealii::DoFHandler objects as other. More...
 
std::shared_ptr< dealii::DoFHandler< dim > > spatial ()
 The underlying spatial dof handler. More...
 
std::shared_ptr< dealii::DoFHandler< 1 > > temporal ()
 The underlying temporal dof handler. More...
 
void distribute_dofs (spacetime::DG_FiniteElement< dim > fe)
 Distribute DoFs in space and time. More...
 
unsigned int n_dofs_spacetime ()
 Total number of space-time degrees of fredom on this slab. More...
 
unsigned int n_dofs_space ()
 Number of spatial degrees of fredom on this slab. More...
 
unsigned int n_dofs_time ()
 Number of temporal degrees of fredom on this slab. More...
 
unsigned int dofs_per_cell_time ()
 Number of temporal dofs in a single element/interval. More...
 
unsigned int dofs_per_cell_space ()
 Number of spatial dofs in a single element. More...
 
spacetime::DG_FiniteElement< dim >::support_type fe_support_type ()
 The underlying support type used for constructing the temporal finite element. More...
 
const dealii::IndexSet & locally_owned_dofs ()
 Return the set of processor local dofs. More...
 

Detailed Description

template<int dim>
class idealii::slab::DoFHandler< dim >

Actual DoFHandler for a specific slab.

This DoFHandler actually handles a spatial and a temporal dealii::DoFHandler object internally.

Currently it is restricted to dG elements in time i.e. spacetime::DG_FiniteElement.

Constructor & Destructor Documentation

◆ DoFHandler() [1/3]

template<int dim>
idealii::slab::DoFHandler< dim >::DoFHandler ( Triangulation< dim > &  tria)

Constructor linking a slab::Triangulation.

Parameters
triaA shared pointer to a slab::Triangulation.

◆ DoFHandler() [2/3]

Constructor linking a parallel::distributed::slab::Triangulation.

Parameters
triaA shared pointer to a parallel::distributed::slab::Triangulation.

◆ DoFHandler() [3/3]

template<int dim>
idealii::slab::DoFHandler< dim >::DoFHandler ( const DoFHandler< dim > &  other)

(shallow) copy constructor. Only the index set and fe support type are actually copied. The underlying pointers will point to the same dealii::DoFHandler objects as other.

Parameters
otherThe DoFHandler to shallow copy.
Warning
This method is mainly needed for adding the DoFHandlers to the spacetime lists and should be used with utmost caution anywhere else.

Member Function Documentation

◆ distribute_dofs()

template<int dim>
void idealii::slab::DoFHandler< dim >::distribute_dofs ( spacetime::DG_FiniteElement< dim >  fe)

Distribute DoFs in space and time.

This function calls the function of the same name of the underlying dof handler objects with the matching spatial or temporal element in fe.

Parameters
feThe spacetime finite element to base the distribution on.

◆ dofs_per_cell_space()

template<int dim>
unsigned int idealii::slab::DoFHandler< dim >::dofs_per_cell_space ( )

Number of spatial dofs in a single element.

Returns
The number of spatial dofs.

◆ dofs_per_cell_time()

template<int dim>
unsigned int idealii::slab::DoFHandler< dim >::dofs_per_cell_time ( )

Number of temporal dofs in a single element/interval.

Returns
The number of temporal dofs i.e. (r+1) for dG(r) elements.

◆ fe_support_type()

template<int dim>
spacetime::DG_FiniteElement<dim>::support_type idealii::slab::DoFHandler< dim >::fe_support_type ( )

The underlying support type used for constructing the temporal finite element.

Returns
The spacetime::DG_FiniteElement<dim>::support_type of the underlying finite element.

◆ locally_owned_dofs()

template<int dim>
const dealii::IndexSet& idealii::slab::DoFHandler< dim >::locally_owned_dofs ( )

Return the set of processor local dofs.

Parallelization is done in space only. Therefore the local space-time dof indices are the local indices of the spatial DoFHandler shifted by the total number of spatial degrees of freedom. Note that the IndexSet is not contiguous and can therefore currently not be used with PETScWrapper classes.

Returns
An IndexSet of the global dof indices owned by the current core.

◆ n_dofs_space()

template<int dim>
unsigned int idealii::slab::DoFHandler< dim >::n_dofs_space ( )

Number of spatial degrees of fredom on this slab.

Returns
The number of dofs based on the spatial finite element and triangulation.

◆ n_dofs_spacetime()

template<int dim>
unsigned int idealii::slab::DoFHandler< dim >::n_dofs_spacetime ( )

Total number of space-time degrees of fredom on this slab.

Returns
The total number of space-time dofs, i.e. n_dofs_space()*n_dofs_time().

◆ n_dofs_time()

template<int dim>
unsigned int idealii::slab::DoFHandler< dim >::n_dofs_time ( )

Number of temporal degrees of fredom on this slab.

Returns
The number of dofs based on the temporal finite element and triangulation.

◆ spatial()

template<int dim>
std::shared_ptr<dealii::DoFHandler<dim> > idealii::slab::DoFHandler< dim >::spatial ( )

The underlying spatial dof handler.

Returns
A shared pointer to the spatial dof handler.

◆ temporal()

template<int dim>
std::shared_ptr<dealii::DoFHandler<1> > idealii::slab::DoFHandler< dim >::temporal ( )

The underlying temporal dof handler.

Returns
A shared pointer to the temporal dof handler.

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