ExaDG
Loading...
Searching...
No Matches
Public Member Functions | Static Public Attributes | List of all members
ExaDG::Acoustics::SpatialOperator< dim, Number > Class Template Reference

#include <spatial_operator.h>

Inheritance diagram for ExaDG::Acoustics::SpatialOperator< dim, Number >:
ExaDG::Acoustics::Interface::SpatialOperator< Number >

Public Member Functions

 SpatialOperator (std::shared_ptr< Grid< dim > const > grid, std::shared_ptr< dealii::Mapping< dim > const > mapping, std::shared_ptr< BoundaryDescriptor< dim > const > boundary_descriptor, std::shared_ptr< FieldFunctions< dim > const > field_functions, Parameters const &parameters, std::string const &field, MPI_Comm const &mpi_comm)
 
void fill_matrix_free_data (MatrixFreeData< dim, Number > &matrix_free_data) const
 
void setup ()
 
void setup (std::shared_ptr< dealii::MatrixFree< dim, Number > const > matrix_free_in, std::shared_ptr< MatrixFreeData< dim, Number > const > matrix_free_data_in)
 
dealii::MatrixFree< dim, Number > const & get_matrix_free () const
 
std::string get_dof_name_pressure () const
 
unsigned int get_dof_index_pressure () const
 
std::string get_dof_name_velocity () const
 
unsigned int get_dof_index_velocity () const
 
unsigned int get_quad_index_pressure () const
 
unsigned int get_quad_index_velocity () const
 
unsigned int get_quad_index_pressure_velocity () const
 
std::shared_ptr< dealii::Mapping< dim > const > get_mapping () const
 
dealii::FiniteElement< dim > const & get_fe_p () const
 
dealii::FiniteElement< dim > const & get_fe_u () const
 
dealii::DoFHandler< dim > const & get_dof_handler_p () const
 
dealii::DoFHandler< dim > const & get_dof_handler_u () const
 
dealii::AffineConstraints< Number > const & get_constraint_p () const
 
dealii::AffineConstraints< Number > const & get_constraint_u () const
 
dealii::types::global_dof_index get_number_of_dofs () const
 
void initialize_dof_vector (BlockVectorType &dst) const final
 
void initialize_dof_vector_pressure (VectorType &dst) const
 
void prescribe_initial_conditions (BlockVectorType &dst, double const time) const final
 
void set_aero_acoustic_source_term (VectorType const &aero_acoustic_source_term_in)
 
void evaluate (BlockVectorType &dst, BlockVectorType const &src, double const time) const final
 
void evaluate_acoustic_operator (BlockVectorType &dst, BlockVectorType const &src, double const time) const
 
void apply_scaled_inverse_mass_operator (BlockVectorType &dst, BlockVectorType const &src) const
 
double calculate_time_step_cfl () const final
 

Static Public Attributes

static unsigned int const block_index_pressure = 0
 
static unsigned int const block_index_velocity = 1
 

Additional Inherited Members

- Public Types inherited from ExaDG::Acoustics::Interface::SpatialOperator< Number >
using BlockVectorType = dealii::LinearAlgebra::distributed::BlockVector<Number>
 

Detailed Description

template<int dim, typename Number>
class ExaDG::Acoustics::SpatialOperator< dim, Number >

Spatial operator to solve for the acoustic conservation equations:

1/c^2 * dp/dt + div (rho * u) = f d(rho * u)/dt + grad p = 0

This operator solves for the pressure p as well as the velocity that is already scaled by the density rho * u.

Member Function Documentation

◆ apply_scaled_inverse_mass_operator()

template<int dim, typename Number >
void ExaDG::Acoustics::SpatialOperator< dim, Number >::apply_scaled_inverse_mass_operator ( BlockVectorType & dst,
BlockVectorType const & src ) const

This function applies the inverse mass matrix and scales the pressure by c^2. This is because the PDE we are solving for reads as

1/c^2 dp/dt + div u = f du/dt + grad p = 0

◆ calculate_time_step_cfl()

template<int dim, typename Number >
double ExaDG::Acoustics::SpatialOperator< dim, Number >::calculate_time_step_cfl ( ) const
finalvirtual

◆ evaluate()

template<int dim, typename Number >
void ExaDG::Acoustics::SpatialOperator< dim, Number >::evaluate ( BlockVectorType & dst,
BlockVectorType const & src,
double const time ) const
finalvirtual

◆ initialize_dof_vector()

template<int dim, typename Number >
void ExaDG::Acoustics::SpatialOperator< dim, Number >::initialize_dof_vector ( BlockVectorType & dst) const
finalvirtual

◆ prescribe_initial_conditions()

template<int dim, typename Number >
void ExaDG::Acoustics::SpatialOperator< dim, Number >::prescribe_initial_conditions ( BlockVectorType & dst,
double const time ) const
finalvirtual

◆ setup() [1/2]

template<int dim, typename Number >
void ExaDG::Acoustics::SpatialOperator< dim, Number >::setup ( )

Call this setup() function if the dealii::MatrixFree object can be set up by the present class.

◆ setup() [2/2]

template<int dim, typename Number >
void ExaDG::Acoustics::SpatialOperator< dim, Number >::setup ( std::shared_ptr< dealii::MatrixFree< dim, Number > const > matrix_free_in,
std::shared_ptr< MatrixFreeData< dim, Number > const > matrix_free_data_in )

Call this setup() function if the dealii::MatrixFree object needs to be created outside this class. The typical use case would be multiphysics-coupling with one MatrixFree object handed over to several single-field solvers. Another typical use case is the use of an ALE formulation.


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