ryujin 2.1.1 revision feb53359f0c9a08baf43c3dfe847d8a9f7d6893a
Namespaces | Classes | Typedefs | Enumerations | Functions | Variables
ryujin Namespace Reference

Namespaces

namespace  Checkpointing
 
namespace  DoFRenumbering
 
namespace  DoFTools
 
namespace  EquationOfStateLibrary
 
namespace  Euler
 
namespace  EulerAEOS
 
namespace  EulerInitialStates
 
namespace  FluxLibrary
 
namespace  Geometries
 
namespace  GridGenerator
 
namespace  internal
 
namespace  Manifolds
 
namespace  NavierStokes
 
namespace  ScalarConservation
 
namespace  ShallowWater
 
namespace  ShallowWaterInitialStates
 
namespace  Skeleton
 

Classes

class  AssemblyCopyData
 
class  AssemblyScratchData
 
class  CubicSpline
 
class  DiagonalPreconditioner
 
class  Discretization
 
class  EquationDispatch
 
class  Geometry
 
struct  get_value_type
 
struct  get_value_type< dealii::VectorizedArray< T, width > >
 
class  HyperbolicModule
 
class  InitialState
 
class  InitialStateLibrary
 
class  InitialStateLibrary< Description, dim, Number >
 
class  InitialValues
 
class  Lazy
 
class  MultiComponentVector
 
class  OfflineData
 
class  ParabolicModule
 
class  Postprocessor
 
class  Quantities
 
class  Restart
 
class  Scope
 
class  SolutionTransfer
 
class  SparseMatrixSIMD
 
class  SparsityPatternSIMD
 
class  StubSolver
 
class  SynchronizationDispatch
 
class  TimeIntegrator
 
class  TimeLoop
 
class  TransfiniteInterpolationManifold
 
class  VTUOutput
 

Typedefs

using Description = Euler::Description
 

Enumerations

enum  Boundary : dealii::types::boundary_id {
  do_nothing = 0 , periodic = 1 , slip = 2 , no_slip = 3 ,
  dirichlet = 4 , dynamic = 5 , dirichlet_momentum = 6
}
 
enum class  Equation {
  euler , euler_aeos , navier_stokes , scalar_conservation ,
  shallow_water
}
 
enum class  IDViolationStrategy { IDViolationStrategy::warn , IDViolationStrategy::raise_exception }
 
enum class  CFLRecoveryStrategy { none , bang_bang_control }
 
enum class  TimeSteppingScheme {
  ssprk_22 , ssprk_33 , erk_11 , erk_22 ,
  erk_33 , erk_43 , erk_54 , strang_ssprk_33_cn ,
  strang_erk_33_cn , strang_erk_43_cn
}
 

Functions

template<int dim, typename Number , typename Callable >
ToFunction< dim, Number, Callable > to_function (const Callable &callable, const unsigned int k)
 
template<typename FT , int problem_dim = FT::dimension, typename TT = typename FT::value_type, typename T = typename TT::value_type>
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, problem_dim, T > contract (const FT &flux_ij, const TT &c_ij)
 
template<typename FT , int problem_dim = FT::dimension>
DEAL_II_ALWAYS_INLINE FT add (const FT &flux_left_ij, const FT &flux_right_ij)
 
void create_parameter_templates (const std::string &parameter_file, const MPI_Comm &mpi_communicator)
 
template<typename Number >
void transform_to_local_range (const dealii::Utilities::MPI::Partitioner &partitioner, dealii::AffineConstraints< Number > &affine_constraints)
 
template<typename VECTOR >
void transform_to_local_range (const dealii::Utilities::MPI::Partitioner &partitioner, VECTOR &vector)
 
std::shared_ptr< const dealii::Utilities::MPI::Partitioner > create_vector_partitioner (const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > &scalar_partitioner, const unsigned int n_components)
 
template<>
float pow (const float x, const float b)
 
template<>
double pow (const double x, const double b)
 
template<>
float fast_pow (const float x, const float b, const Bias)
 
template<>
double fast_pow (const double x, const double b, const Bias)
 
template<typename T >
fast_pow_impl (const T x, const T b, const Bias)
 
void print_revision_and_version (std::ostream &stream)
 
void print_compile_time_parameters (std::ostream &stream)
 
Functions
template<typename Number >
DEAL_II_ALWAYS_INLINE void quadratic_newton_step (Number &p_1, Number &p_2, const Number phi_p_1, const Number phi_p_2, const Number dphi_p_1, const Number dphi_p_2, const Number sign=Number(1.0))
 
SIMD based access to vectors and arrays of vectors
template<typename T , typename V >
DEAL_II_ALWAYS_INLINE T get_entry (const V &vector, unsigned int i)
 
template<typename T , typename T2 >
DEAL_II_ALWAYS_INLINE T get_entry (const std::vector< T2 > &vector, unsigned int i)
 
template<typename T , typename V >
DEAL_II_ALWAYS_INLINE T get_entry (const V &vector, const unsigned int *js)
 
template<typename T , typename T2 >
DEAL_II_ALWAYS_INLINE T get_entry (const std::vector< T2 > &vector, const unsigned int *js)
 
template<typename T , typename V >
DEAL_II_ALWAYS_INLINE void write_entry (V &vector, const T &values, unsigned int i)
 
template<typename T , typename T2 >
DEAL_II_ALWAYS_INLINE void write_entry (std::vector< T2 > &vector, const T &values, unsigned int i)
 
template<int rank, int dim, std::size_t width, typename Number >
DEAL_II_ALWAYS_INLINE dealii::Tensor< rank, dim, Number > serialize_tensor (const dealii::Tensor< rank, dim, dealii::VectorizedArray< Number, width > > &vectorized, const unsigned int k)
 
template<int rank, int dim, typename Number >
DEAL_II_ALWAYS_INLINE dealii::Tensor< rank, dim, Number > serialize_tensor (const dealii::Tensor< rank, dim, Number > &serial, const unsigned int k)
 
template<int rank, int dim, std::size_t width, typename Number >
DEAL_II_ALWAYS_INLINE void assign_serial_tensor (dealii::Tensor< rank, dim, dealii::VectorizedArray< Number, width > > &result, const dealii::Tensor< rank, dim, Number > &serial, const unsigned int k)
 
template<int rank, int dim, typename Number >
DEAL_II_ALWAYS_INLINE void assign_serial_tensor (dealii::Tensor< rank, dim, Number > &result, const dealii::Tensor< rank, dim, Number > &serial, const unsigned int k)
 

Variables

template<int dim>
constexpr bool have_distributed_triangulation
 

Transcendental and other mathematical operations

enum class  Bias { none , max , min }
 
template<typename Number >
DEAL_II_ALWAYS_INLINE Number positive_part (const Number number)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE Number negative_part (const Number number)
 
template<int N, typename T >
fixed_power (const T x)
 
template<typename T >
pow (const T x, const T b)
 
template<typename T , std::size_t width>
dealii::VectorizedArray< T, width > pow (const dealii::VectorizedArray< T, width > x, const T b)
 
template<typename T , std::size_t width>
dealii::VectorizedArray< T, width > pow (const dealii::VectorizedArray< T, width > x, const dealii::VectorizedArray< T, width > b)
 
template<typename T >
fast_pow (const T x, const T b, const Bias bias=Bias::none)
 
template<typename T , std::size_t width>
dealii::VectorizedArray< T, width > fast_pow (const dealii::VectorizedArray< T, width > x, const T b, const Bias bias=Bias::none)
 
template<typename T , std::size_t width>
dealii::VectorizedArray< T, width > fast_pow (const dealii::VectorizedArray< T, width > x, const dealii::VectorizedArray< T, width > b, const Bias bias=Bias::none)
 

Type traits and packed index handling

template<typename T >
unsigned int get_stride_size = 1
 
template<unsigned int length, typename Functor >
DEAL_II_ALWAYS_INLINE auto generate_iterators (Functor f) -> std::array< auto, length >
 
template<typename T >
DEAL_II_ALWAYS_INLINE void increment_iterators (T &iterators)
 

Typedef Documentation

◆ Description

Definition at line 15 of file initial_state_library.template.h.

Enumeration Type Documentation

◆ Equation

enum class ryujin::Equation
strong

An enum class that controls what problem description to use.

Enumerator
euler 

The compressible Euler equations of gas dynamics with a polytropic gas equation of state

euler_aeos 

The compressible Euler equations of gas dynamics with arbitrary equation of state

navier_stokes 

The compressible Navier-Stokes equations of gas dynamics with a polytropic gas equation of state, Newtonian fluid viscosity model, and a heat flux governed by Fourier's law.

scalar_conservation 

A scalar conservation equation with a user-specified flux depending on the state.

shallow_water 

The shallow water equations

Definition at line 26 of file equation_dispatch.h.

◆ Bias

enum class ryujin::Bias
strong

Controls the bias of the fast_pow() functions.

Enumerator
none 

No specific bias.

max 

Guarantee an upper bound, i.e., fast_pow(x,b) >= pow(x,b) provided that FIXME

min 

Guarantee a lower bound, i.e., fast_pow(x,b) >= pow(x,b) provided that FIXME

Definition at line 177 of file simd.h.

◆ CFLRecoveryStrategy

enum class ryujin::CFLRecoveryStrategy
strong

Controls the chosen invariant domain / CFL recovery strategy.

Enumerator
none 

Step with the chosen "cfl max" value and do nothing in case an invariant domain and or CFL condition violation is detected.

bang_bang_control 

Step with the chosen "cfl max" value and, in case an invariant domain and or CFL condition violation is detected, the time step is repeated with "cfl min". If this is unsuccessful as well, a warning is emitted.

Definition at line 21 of file time_integrator.h.

◆ TimeSteppingScheme

enum class ryujin::TimeSteppingScheme
strong

Controls the chosen time-stepping scheme.

Enumerator
ssprk_22 

The strong stability preserving Runge Kutta method of order 2, SSPRK(2,2;1/2), with the following butcher tableau

\begin{align*} \begin{array}{c|ccc} 0 & 0 \\ \tfrac{1}{2} & \tfrac{1}{2} & 0 \\ \hline 1 & 1 & 0 \end{array} \end{align*}

ssprk_33 

The strong stability preserving Runge Kutta method of order 3, SSPRK(3,3;1/3), with the following butcher tableau

\begin{align*} \begin{array}{c|ccc} 0 & 0 \\ 1 & 1 & 0 \\ \tfrac{1}{2} & \tfrac{1}{4} & \tfrac{1}{4} & 0\\ \hline 1 & \tfrac{1}{6} & \tfrac{1}{6} & \tfrac{2}{3} \end{array} \end{align*}

erk_11 

The explicit Runge-Kutta method RK(1,1;1), aka a simple, forward Euler step.

erk_22 

The explicit Runge-Kutta method RK(2,2;1) with the butcher tableau

\begin{align*} \begin{array}{c|ccc} 0 & 0 \\ \tfrac{1}{2} & \tfrac{1}{2} & 0 \\ \hline 1 & 0 & 1 \end{array} \end{align*}

erk_33 

The explicit Runge-Kutta method RK(3,3;1) with the butcher tableau

\begin{align*} \begin{array}{c|ccc} 0 & 0 \\ \tfrac{1}{3} & \tfrac{1}{3} & 0 \\ \tfrac{2}{3} & 0 & \tfrac{2}{3} & 0 \\ \hline 1 & \tfrac{1}{4} & 0 & \tfrac{3}{4} \end{array} \end{align*}

erk_43 

The explicit Runge-Kutta method RK(4,3;1) with the butcher tableau

\begin{align*} \begin{array}{c|ccc} 0 & 0 \\ \tfrac{1}{4} & \tfrac{1}{4} & 0 \\ \tfrac{1}{2} & 0 & \tfrac{1}{2} & 0 \\ \tfrac{3}{4} & 0 & \tfrac{1}{4} & \tfrac{1}{2} & 0 \\ \hline 1 & 0 & \tfrac{2}{3} & -\tfrac{1}{3} & \tfrac{2}{3} \end{array} \end{align*}

erk_54 

The explicit Runge-Kutta method RK(5,4;1) with the butcher tableau TODO

strang_ssprk_33_cn 

A Strang split using ssprk 33 for the hyperbolic subproblem and Crank-Nicolson for the parabolic subproblem

strang_erk_33_cn 

A Strang split using erk 33 for the hyperbolic subproblem and Crank-Nicolson for the parabolic subproblem

strang_erk_43_cn 

A Strang split using erk 43 for the hyperbolic subproblem and Crank-Nicolson for the parabolic subproblem

Definition at line 41 of file time_integrator.h.

Function Documentation

◆ contract()

template<typename FT , int problem_dim = FT::dimension, typename TT = typename FT::value_type, typename T = typename TT::value_type>
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, problem_dim, T > ryujin::contract ( const FT &  flux_ij,
const TT &  c_ij 
)
inline

◆ add()

template<typename FT , int problem_dim = FT::dimension>
DEAL_II_ALWAYS_INLINE FT ryujin::add ( const FT &  flux_left_ij,
const FT &  flux_right_ij 
)
inline

◆ create_parameter_templates()

void ryujin::create_parameter_templates ( const std::string &  parameter_file,
const MPI_Comm &  mpi_communicator 
)

Definition at line 241 of file equation_dispatch.h.

◆ pow() [1/2]

template<>
float ryujin::pow ( const float  x,
const float  b 
)

Definition at line 218 of file simd.template.h.

References pow().

◆ pow() [2/2]

template<>
double ryujin::pow ( const double  x,
const double  b 
)

Definition at line 227 of file simd.template.h.

References pow().

Referenced by fast_pow(), and pow().

◆ fast_pow() [1/2]

template<>
float ryujin::fast_pow ( const float  x,
const float  b,
const  Bias 
)

Definition at line 280 of file simd.template.h.

References pow().

◆ fast_pow() [2/2]

template<>
double ryujin::fast_pow ( const double  x,
const double  b,
const  Bias 
)

Definition at line 289 of file simd.template.h.

References pow().

◆ fast_pow_impl()

template<typename T >
T ryujin::fast_pow_impl ( const T  x,
const T  b,
const  Bias 
)

Definition at line 25 of file simd_fast_pow.template.h.

References pow().

Referenced by fast_pow().

◆ print_revision_and_version()

void ryujin::print_revision_and_version ( std::ostream &  stream)

Print git revision and version info.

Referenced by ryujin::TimeLoop< Description, dim, Number >::print_parameters().

◆ print_compile_time_parameters()

void ryujin::print_compile_time_parameters ( std::ostream &  stream)

Print compile time parameters.

Referenced by ryujin::TimeLoop< Description, dim, Number >::print_parameters().

Variable Documentation

◆ have_distributed_triangulation

template<int dim>
constexpr bool ryujin::have_distributed_triangulation
constexpr
Initial value:
=
std::is_same<typename Discretization<dim>::Triangulation,
dealii::parallel::distributed::Triangulation<dim>>::value

A templated constexpr boolean that is true if we use a parallel distributed triangulation (for the specified dimension).

Definition at line 156 of file discretization.h.