[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
PHiLiP::Physics::PhysicsBase< dim, nstate, real > Class Template Referenceabstract

Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived. More...

#include <physics.h>

Inheritance diagram for PHiLiP::Physics::PhysicsBase< dim, nstate, real >:
Collaboration diagram for PHiLiP::Physics::PhysicsBase< dim, nstate, real >:

Public Types

using NonPhysicalBehaviorEnum = Parameters::AllParameters::NonPhysicalBehaviorEnum
 

Public Member Functions

 PhysicsBase (const Parameters::AllParameters *const parameters_input, const bool has_nonzero_diffusion_input, const bool has_nonzero_physical_source_input, const dealii::Tensor< 2, 3, double > input_diffusion_tensor=Parameters::ManufacturedSolutionParam::get_default_diffusion_tensor(), std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function_input=nullptr)
 Default constructor that will set the constants.
 
 PhysicsBase (const Parameters::AllParameters *const parameters_input, const bool has_nonzero_diffusion_input, const bool has_nonzero_physical_source_input, std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function_input=nullptr)
 Constructor that will call default constructor.
 
virtual ~PhysicsBase ()=default
 Virtual destructor required for abstract classes.
 
virtual std::array< dealii::Tensor< 1, dim, real >, nstate > convective_flux (const std::array< real, nstate > &solution) const =0
 Convective fluxes that will be differentiated once in space.
 
virtual std::array< dealii::Tensor< 1, dim, real >, nstate > convective_numerical_split_flux (const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &conservative_soln2) const
 Convective Numerical Split Flux for split form.
 
virtual std::array< real, nstate > compute_entropy_variables (const std::array< real, nstate > &conservative_soln) const =0
 Computes the entropy variables.
 
virtual std::array< real, nstate > compute_conservative_variables_from_entropy_variables (const std::array< real, nstate > &entropy_var) const =0
 Computes the conservative variables from the entropy variables.
 
virtual std::array< real, nstate > convective_eigenvalues (const std::array< real, nstate > &, const dealii::Tensor< 1, dim, real > &) const =0
 Spectral radius of convective term Jacobian. More...
 
virtual real max_convective_eigenvalue (const std::array< real, nstate > &soln) const =0
 Maximum convective eigenvalue.
 
virtual real max_convective_normal_eigenvalue (const std::array< real, nstate > &soln, const dealii::Tensor< 1, dim, real > &normal) const
 Maximum convective normal eigenvalue (used in Lax-Friedrichs)
 
virtual real max_viscous_eigenvalue (const std::array< real, nstate > &soln) const =0
 Maximum viscous eigenvalue.
 
virtual std::array< dealii::Tensor< 1, dim, real >, nstate > dissipative_flux (const std::array< real, nstate > &solution, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const =0
 Dissipative fluxes that will be differentiated ONCE in space.
 
virtual std::array< real, nstate > source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &solution, const real current_time, const dealii::types::global_dof_index cell_index) const =0
 Artificial dissipative fluxes that will be differentiated ONCE in space. More...
 
virtual std::array< real, nstate > physical_source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &solution, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const
 Physical source term that does require differentiation.
 
virtual std::array< real, nstate > artificial_source_term (const real viscosity_coefficient, const dealii::Point< dim, real > &pos, const std::array< real, nstate > &solution) const
 Artificial source term that does not require differentiation stemming from artificial dissipation.
 
virtual void boundary_face_values (const int, const dealii::Point< dim, real > &, const dealii::Tensor< 1, dim, real > &, const std::array< real, nstate > &, const std::array< dealii::Tensor< 1, dim, real >, nstate > &, std::array< real, nstate > &, std::array< dealii::Tensor< 1, dim, real >, nstate > &) const =0
 Evaluates boundary values and gradients on the other side of the face.
 
virtual dealii::Vector< double > post_compute_derived_quantities_vector (const dealii::Vector< double > &uh, const std::vector< dealii::Tensor< 1, dim > > &, const std::vector< dealii::Tensor< 2, dim > > &, const dealii::Tensor< 1, dim > &, const dealii::Point< dim > &) const
 Returns current vector solution to be used by PhysicsPostprocessor to output current solution. More...
 
virtual dealii::Vector< double > post_compute_derived_quantities_scalar (const double &uh, const dealii::Tensor< 1, dim > &, const dealii::Tensor< 2, dim > &, const dealii::Tensor< 1, dim > &, const dealii::Point< dim > &) const
 Returns current scalar solution to be used by PhysicsPostprocessor to output current solution. More...
 
virtual std::vector< std::string > post_get_names () const
 Returns names of the solution to be used by PhysicsPostprocessor to output current solution. More...
 
virtual std::vector< dealii::DataComponentInterpretation::DataComponentInterpretation > post_get_data_component_interpretation () const
 Returns DataComponentInterpretation of the solution to be used by PhysicsPostprocessor to output current solution. More...
 
virtual dealii::UpdateFlags post_get_needed_update_flags () const
 Returns required update flags of the solution to be used by PhysicsPostprocessor to output current solution. More...
 
template<typename real2 >
real2 handle_non_physical_result (const std::string message="") const
 Function to handle nonphysical results. More...
 

Public Attributes

const bool has_nonzero_diffusion
 Flag to signal that diffusion term is non-zero.
 
const bool has_nonzero_physical_source
 Flag to signal that physical source term is non-zero.
 
const Parameters::AllParameters *const all_parameters
 Pointer to parameters object.
 
const NonPhysicalBehaviorEnum non_physical_behavior_type
 Determines type of nonphysical behavior.
 
std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function
 Manufactured solution function.
 
const double BIG_NUMBER = 1e100
 BIG_NUMBER which is returned in place of NaN according to handle_non_physical_result() More...
 

Protected Attributes

dealii::ConditionalOStream pcout
 ConditionalOStream. More...
 
dealii::Tensor< 2, dim, double > diffusion_tensor
 Anisotropic diffusion matrix. More...
 

Detailed Description

template<int dim, int nstate, typename real>
class PHiLiP::Physics::PhysicsBase< dim, nstate, real >

Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived.

Main interface for all the convective and diffusive terms.

LinearAdvection, Diffusion, ConvectionDiffusion, Euler are derived from this class.

Partial differential equation is given by the divergence of the convective and diffusive flux equal to the source term

\[ \boldsymbol{\nabla} \cdot ( \mathbf{F}_{conv}( u ) + \mathbf{F}_{diss}( u, \boldsymbol{\nabla}(u) ) = s(\mathbf{x}) \]

Definition at line 34 of file physics.h.

Member Function Documentation

◆ convective_eigenvalues()

◆ handle_non_physical_result()

template<int dim, int nstate, typename real >
template<typename real2 >
template FadType PHiLiP::Physics::PhysicsBase< dim, nstate, real >::handle_non_physical_result< FadType > ( const std::string  message = "") const

Function to handle nonphysical results.

This is to be called by derived physics classes when a nonphysical quantity is detected there. e.g., negative density in Euler. The behavior is determined by the value of non_physical_behavior_type. Returns BIG_NUMBER.

Definition at line 208 of file physics.cpp.

◆ post_compute_derived_quantities_scalar()

template<int dim, int nstate, typename real >
dealii::Vector< double > PHiLiP::Physics::PhysicsBase< dim, nstate, real >::post_compute_derived_quantities_scalar ( const double &  uh,
const dealii::Tensor< 1, dim > &  ,
const dealii::Tensor< 2, dim > &  ,
const dealii::Tensor< 1, dim > &  ,
const dealii::Point< dim > &   
) const
virtual

Returns current scalar solution to be used by PhysicsPostprocessor to output current solution.

The implementation in this Physics base class simply returns the stored solution.

Definition at line 160 of file physics.cpp.

◆ post_compute_derived_quantities_vector()

template<int dim, int nstate, typename real >
dealii::Vector< double > PHiLiP::Physics::PhysicsBase< dim, nstate, real >::post_compute_derived_quantities_vector ( const dealii::Vector< double > &  uh,
const std::vector< dealii::Tensor< 1, dim > > &  ,
const std::vector< dealii::Tensor< 2, dim > > &  ,
const dealii::Tensor< 1, dim > &  ,
const dealii::Point< dim > &   
) const
virtual

◆ post_get_data_component_interpretation()

template<int dim, int nstate, typename real >
std::vector< dealii::DataComponentInterpretation::DataComponentInterpretation > PHiLiP::Physics::PhysicsBase< dim, nstate, real >::post_get_data_component_interpretation ( ) const
virtual

◆ post_get_names()

template<int dim, int nstate, typename real >
std::vector< std::string > PHiLiP::Physics::PhysicsBase< dim, nstate, real >::post_get_names ( ) const
virtual

◆ post_get_needed_update_flags()

template<int dim, int nstate, typename real >
dealii::UpdateFlags PHiLiP::Physics::PhysicsBase< dim, nstate, real >::post_get_needed_update_flags ( ) const
virtual

◆ source_term()

template<int dim, int nstate, typename real>
virtual std::array<real,nstate> PHiLiP::Physics::PhysicsBase< dim, nstate, real >::source_term ( const dealii::Point< dim, real > &  pos,
const std::array< real, nstate > &  solution,
const real  current_time,
const dealii::types::global_dof_index  cell_index 
) const
pure virtual

Member Data Documentation

◆ BIG_NUMBER

template<int dim, int nstate, typename real>
const double PHiLiP::Physics::PhysicsBase< dim, nstate, real >::BIG_NUMBER = 1e100

BIG_NUMBER which is returned in place of NaN according to handle_non_physical_result()

Type double so that typecasting works with all real types

Definition at line 206 of file physics.h.

◆ diffusion_tensor

template<int dim, int nstate, typename real>
dealii::Tensor<2,dim,double> PHiLiP::Physics::PhysicsBase< dim, nstate, real >::diffusion_tensor
protected

Anisotropic diffusion matrix.

As long as the diagonal components are positive and diagonally dominant we should have a stable diffusive system

Definition at line 218 of file physics.h.

◆ pcout

template<int dim, int nstate, typename real>
dealii::ConditionalOStream PHiLiP::Physics::PhysicsBase< dim, nstate, real >::pcout
protected

ConditionalOStream.

Used as std::cout, but only prints if mpi_rank == 0

Definition at line 212 of file physics.h.


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