ryujin 2.1.1 revision 61395bc58278f41c93a581f9251f0a9c9dc4b8a9
Public Types | Public Member Functions | List of all members
ryujin::OfflineData< dim, Number > Class Template Reference

#include <source/offline_data.h>

Inheritance diagram for ryujin::OfflineData< dim, Number >:
Inheritance graph
[legend]
Collaboration diagram for ryujin::OfflineData< dim, Number >:
Collaboration graph
[legend]

Public Types

using scalar_type = dealii::LinearAlgebra::distributed::Vector< Number >
 
using boundary_description = std::tuple< dealii::Tensor< 1, dim, Number >, Number, Number, dealii::types::boundary_id, dealii::Point< dim > >
 

Public Member Functions

 OfflineData (const MPI_Comm &mpi_communicator, const Discretization< dim > &discretization, const std::string &subsection="/OfflineData")
 
void prepare (const unsigned int problem_dimension)
 
auto & dof_handler () const
 
auto & affine_constraints () const
 
const auto & scalar_partitioner () const
 
const auto & vector_partitioner () const
 
auto & n_export_indices () const
 
auto & n_locally_internal () const
 
auto & n_locally_owned () const
 
auto & n_locally_relevant () const
 
auto & boundary_map () const
 
auto & coupling_boundary_pairs () const
 
auto & level_boundary_map () const
 
auto & sparsity_pattern () const
 
auto & sparsity_pattern_simd () const
 
auto & mass_matrix () const
 
auto & lumped_mass_matrix () const
 
auto & lumped_mass_matrix_inverse () const
 
auto & level_lumped_mass_matrix () const
 
auto & betaij_matrix () const
 
auto & cij_matrix () const
 
auto & measure_of_omega () const
 
auto & discretization () const
 
template<typename ITERATOR1 , typename ITERATOR2 >
OfflineData< dim, Number >::boundary_map_type construct_boundary_map (const ITERATOR1 &begin, const ITERATOR2 &end, const Utilities::MPI::Partitioner &partitioner) const
 
template<typename ITERATOR1 , typename ITERATOR2 >
OfflineData< dim, Number >::coupling_boundary_pairs_type collect_coupling_boundary_pairs (const ITERATOR1 &begin, const ITERATOR2 &end, const Utilities::MPI::Partitioner &partitioner) const
 

Detailed Description

template<int dim, typename Number = double>
class ryujin::OfflineData< dim, Number >

A class to store all data that can be precomputed offline.

This class takes a reference to a Discretization object (that itself holds a Triangulation, FiniteElement, Mapping, and Quadrature object).

Most notably this class sets up a DoFHandler, the SparsityPattern, various IndexSet objects to hold locally owned and locally relevant indices, and precomputes all matrices (mass matrix, lumped mass matrix, $c_{ij}$ matrices, and $n_{ij}$ matrices).

After prepare() is called, all getter functions return valid references.

Note
The offline data precomputed in this class is problem independent, it only depends on the chosen geometry and ansatz stored in the Discretization class.

Definition at line 47 of file offline_data.h.

Member Typedef Documentation

◆ scalar_type

template<int dim, typename Number = double>
using ryujin::OfflineData< dim, Number >::scalar_type = dealii::LinearAlgebra::distributed::Vector<Number>

Shorthand typedef for dealii::LinearAlgebra::distributed::Vector<Number>.

Definition at line 54 of file offline_data.h.

◆ boundary_description

template<int dim, typename Number = double>
using ryujin::OfflineData< dim, Number >::boundary_description = std::tuple<dealii::Tensor<1, dim, Number> , Number , Number , dealii::types::boundary_id , dealii::Point<dim> >

A tuple describing global dof index, boundary normal, normal mass, boundary mass, boundary id, and position of the boundary degree of freedom.

Definition at line 61 of file offline_data.h.

Constructor & Destructor Documentation

◆ OfflineData()

template<int dim, typename Number >
ryujin::OfflineData< dim, Number >::OfflineData ( const MPI_Comm &  mpi_communicator,
const Discretization< dim > &  discretization,
const std::string &  subsection = "/OfflineData< dim, Number >" 
)

Constructor

Definition at line 37 of file offline_data.template.h.

Member Function Documentation

◆ prepare()

template<int dim, typename Number = double>
void ryujin::OfflineData< dim, Number >::prepare ( const unsigned int  problem_dimension)
inline

Prepare offline data. A call to prepare() internally calls setup() and assemble().

The problem_dimension parameter is used to setup up an appropriately sized vector partitioner for the MultiComponentVector.

Definition at line 82 of file offline_data.h.

◆ dof_handler()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::dof_handler ( ) const
inline

The DofHandler for our (scalar) CG ansatz space in (deal.II typical) global numbering.

Definition at line 93 of file offline_data.h.

Referenced by ryujin::Checkpointing::load_state_vector(), and ryujin::Checkpointing::write_checkpoint().

◆ affine_constraints()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::affine_constraints ( ) const
inline

An AffineConstraints object storing constraints in (Deal.II typical) global numbering.

Definition at line 99 of file offline_data.h.

◆ scalar_partitioner()

template<int dim, typename Number = double>
const auto & ryujin::OfflineData< dim, Number >::scalar_partitioner ( ) const
inline

An MPI partitioner for all parallel distributed vectors storing a scalar quantity.

Definition at line 105 of file offline_data.h.

Referenced by ryujin::Checkpointing::load_state_vector(), and ryujin::Checkpointing::write_checkpoint().

◆ vector_partitioner()

template<int dim, typename Number = double>
const auto & ryujin::OfflineData< dim, Number >::vector_partitioner ( ) const
inline

An MPI partitioner for all parallel distributed vectors storing a vector-valued quantity of size HyperbolicSystem::problem_dimension.

Definition at line 112 of file offline_data.h.

◆ n_export_indices()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::n_export_indices ( ) const
inline

The subinterval \([0,\texttt{n_export_indices()})\) contains all (SIMD-vectorized) indices of the interval \([0,\texttt{n_locally_internal()})\) that are exported to neighboring MPI ranks.

Note
The interval \([\texttt{n_locally_internal()}, \texttt{n_locally_relevant()})\) (consisting of non-SIMD-vectorized indices) contains additional degrees of freedom that might have to be exported to neighboring MPI ranks.

Definition at line 125 of file offline_data.h.

◆ n_locally_internal()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::n_locally_internal ( ) const
inline

Number of locally owned internal degrees of freedom: In (MPI rank) local numbering all indices in the half open interval [0, n_locally_internal_) are owned by this processor, have standard connectivity, and are not situated at a boundary.

Definition at line 133 of file offline_data.h.

◆ n_locally_owned()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::n_locally_owned ( ) const
inline

Number of locally owned degrees of freedom: In (MPI rank) local numbering all indices in the half open interval [0, n_locally_owned_) are owned by this processor.

Definition at line 140 of file offline_data.h.

◆ n_locally_relevant()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::n_locally_relevant ( ) const
inline

Number of locally relevant degrees of freedom: This number is the toal number of degrees of freedom we store locally on this MPI rank. I.e., we can access the half open interval [0, n_locally_relevant_) on this machine.

Definition at line 148 of file offline_data.h.

◆ boundary_map()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::boundary_map ( ) const
inline

The boundary map. Local numbering.

For every degree of freedom that has nonzero support at the boundary we record the global degree of freedom index along with a weighted boundary normal, the associated boundary id, and position.

This map is later used in OfflineData to handle boundary degrees of freedom after every time step (for example to implement reflective boundary conditions).

Definition at line 161 of file offline_data.h.

Referenced by ryujin::NavierStokes::VelocityMatrix< dim, Number, Number2 >::vmult(), and ryujin::NavierStokes::EnergyMatrix< dim, Number, Number2 >::vmult().

◆ coupling_boundary_pairs()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::coupling_boundary_pairs ( ) const
inline

A vector of tuples describing coupling degrees of freedom i and j where both degrees of freedom are collocated at the boundary (and hence the d_ij matrix has to be symmetrized). The function returns a reference to a vector of tuples consisting of (i, col_idx, j).

Definition at line 169 of file offline_data.h.

◆ level_boundary_map()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::level_boundary_map ( ) const
inline

◆ sparsity_pattern()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::sparsity_pattern ( ) const
inline

A sparsity pattern for (standard deal.II) matrices storing indices in (Deal.II typical) global numbering.

Definition at line 181 of file offline_data.h.

◆ sparsity_pattern_simd()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::sparsity_pattern_simd ( ) const
inline

A sparsity pattern for matrices in vectorized format. Local numbering.

Definition at line 187 of file offline_data.h.

◆ mass_matrix()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::mass_matrix ( ) const
inline

The mass matrix. (SIMD storage, local numbering)

Definition at line 192 of file offline_data.h.

◆ lumped_mass_matrix()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::lumped_mass_matrix ( ) const
inline

◆ lumped_mass_matrix_inverse()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::lumped_mass_matrix_inverse ( ) const
inline

The inverse of the lumped mass matrix.

Definition at line 202 of file offline_data.h.

◆ level_lumped_mass_matrix()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::level_lumped_mass_matrix ( ) const
inline

◆ betaij_matrix()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::betaij_matrix ( ) const
inline

The stiffness matrix \((beta_{ij})\): \(\beta_{ij} = \nabla\varphi_{j}\cdot\nabla\varphi_{i}\) (SIMD storage, local numbering)

Definition at line 215 of file offline_data.h.

◆ cij_matrix()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::cij_matrix ( ) const
inline

The \((c_{ij})\) matrix. (SIMD storage, local numbering)

Definition at line 220 of file offline_data.h.

◆ measure_of_omega()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::measure_of_omega ( ) const
inline

Size of computational domain.

Definition at line 225 of file offline_data.h.

◆ discretization()

template<int dim, typename Number = double>
auto & ryujin::OfflineData< dim, Number >::discretization ( ) const
inline

Returns a reference of the underlying Discretization object.

Definition at line 230 of file offline_data.h.

Referenced by ryujin::Checkpointing::write_checkpoint().

◆ construct_boundary_map()

template<int dim, typename Number = double>
template<typename ITERATOR1 , typename ITERATOR2 >
OfflineData< dim, Number >::boundary_map_type ryujin::OfflineData< dim, Number >::construct_boundary_map ( const ITERATOR1 &  begin,
const ITERATOR2 &  end,
const Utilities::MPI::Partitioner &  partitioner 
) const

Definition at line 695 of file offline_data.template.h.

References ryujin::periodic.

◆ collect_coupling_boundary_pairs()

template<int dim, typename Number = double>
template<typename ITERATOR1 , typename ITERATOR2 >
OfflineData< dim, Number >::coupling_boundary_pairs_type ryujin::OfflineData< dim, Number >::collect_coupling_boundary_pairs ( const ITERATOR1 &  begin,
const ITERATOR2 &  end,
const Utilities::MPI::Partitioner &  partitioner 
) const

Definition at line 850 of file offline_data.template.h.

References ryujin::periodic.


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