[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Physics Model equations. Derived from PhysicsBase, holds a baseline physics and model terms and equations. More...
#include <physics_model.h>
Public Member Functions | |
PhysicsModel (const Parameters::AllParameters *const parameters_input, Parameters::AllParameters::PartialDifferentialEquation baseline_physics_type, std::shared_ptr< ModelBase< dim, nstate, real > > model_input, std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function, const bool has_nonzero_diffusion, const bool has_nonzero_physical_source) | |
Constructor. | |
std::array< dealii::Tensor< 1, dim, real >, nstate > | convective_flux (const std::array< real, nstate > &conservative_soln) const |
Convective flux: \( \mathbf{F}_{conv} \). | |
std::array< dealii::Tensor< 1, dim, real >, nstate > | dissipative_flux (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const |
Dissipative (i.e. viscous) flux: \( \mathbf{F}_{diss} \). | |
std::array< real, nstate > | physical_source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const |
Physical source term. | |
std::array< real, nstate > | source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_soln, const real current_time, const dealii::types::global_dof_index cell_index) const |
Source term that does not require differentiation. | |
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. | |
std::array< real, nstate > | compute_entropy_variables (const std::array< real, nstate > &conservative_soln) const |
Computes the entropy variables. | |
std::array< real, nstate > | compute_conservative_variables_from_entropy_variables (const std::array< real, nstate > &entropy_var) const |
Computes the conservative variables from the entropy variables. | |
std::array< real, nstate > | convective_eigenvalues (const std::array< real, nstate > &, const dealii::Tensor< 1, dim, real > &) const |
real | max_convective_eigenvalue (const std::array< real, nstate > &soln) const |
Maximum convective eigenvalue. | |
real | max_convective_normal_eigenvalue (const std::array< real, nstate > &soln, const dealii::Tensor< 1, dim, real > &normal) const override |
Maximum convective normal eigenvalue (used in Lax-Friedrichs) | |
real | max_viscous_eigenvalue (const std::array< real, nstate > &soln) const |
Maximum viscous eigenvalue. | |
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 |
Evaluates boundary values and gradients on the other side of the face. | |
dealii::Vector< double > | post_compute_derived_quantities_vector (const dealii::Vector< double > &uh, const std::vector< dealii::Tensor< 1, dim > > &duh, const std::vector< dealii::Tensor< 2, dim > > &dduh, const dealii::Tensor< 1, dim > &normals, const dealii::Point< dim > &evaluation_points) const |
Returns current vector solution to be used by PhysicsPostprocessor to output current solution. More... | |
std::vector< std::string > | post_get_names () const |
Returns names of the solution to be used by PhysicsPostprocessor to output current solution. More... | |
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... | |
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... | |
![]() | |
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< 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 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... | |
template<typename real2 > | |
real2 | handle_non_physical_result (const std::string message="") const |
Function to handle nonphysical results. More... | |
Public Attributes | |
const int | n_model_equations |
Number of model equations (i.e. those additional to the baseline physics) | |
std::shared_ptr< PhysicsBase< dim, nstate_baseline_physics, real > > | physics_baseline |
Baseline physics object with nstate==nstate_baseline_physics. | |
std::shared_ptr< ModelBase< dim, nstate, real > > | model |
Model object. | |
![]() | |
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 | |
const MPI_Comm | mpi_communicator |
MPI communicator. | |
const int | mpi_rank |
MPI rank. | |
const int | n_mpi |
dealii::ConditionalOStream | pcout |
ConditionalOStream. More... | |
![]() | |
dealii::ConditionalOStream | pcout |
ConditionalOStream. More... | |
dealii::Tensor< 2, dim, double > | diffusion_tensor |
Anisotropic diffusion matrix. More... | |
Additional Inherited Members | |
![]() | |
using | NonPhysicalBehaviorEnum = Parameters::AllParameters::NonPhysicalBehaviorEnum |
Physics Model equations. Derived from PhysicsBase, holds a baseline physics and model terms and equations.
Definition at line 14 of file physics_model.h.
|
virtual |
Spectral radius of convective term Jacobian. Used for scalar dissipation
Implements PHiLiP::Physics::PhysicsBase< dim, nstate, real >.
Definition at line 219 of file physics_model.cpp.
|
virtual |
Returns current vector solution to be used by PhysicsPostprocessor to output current solution.
The implementation in this Physics base class simply returns the stored solution.
Reimplemented from PHiLiP::Physics::PhysicsBase< dim, nstate, real >.
Definition at line 341 of file physics_model.cpp.
|
virtual |
Returns DataComponentInterpretation of the solution to be used by PhysicsPostprocessor to output current solution.
Treats every solution state as an independent scalar.
Reimplemented from PHiLiP::Physics::PhysicsBase< dim, nstate, real >.
Definition at line 390 of file physics_model.cpp.
|
virtual |
Returns names of the solution to be used by PhysicsPostprocessor to output current solution.
The implementation in this Physics base class simply returns "state0, state1, etc.".
Reimplemented from PHiLiP::Physics::PhysicsBase< dim, nstate, real >.
Definition at line 374 of file physics_model.cpp.
|
virtual |
Returns required update flags of the solution to be used by PhysicsPostprocessor to output current solution.
Only update the solution at the output points.
Reimplemented from PHiLiP::Physics::PhysicsBase< dim, nstate, real >.
Definition at line 407 of file physics_model.cpp.
|
protected |
Number of MPI processes.
Definition at line 131 of file physics_model.h.
|
protected |
ConditionalOStream.
Used as std::cout, but only prints if mpi_rank == 0
Definition at line 135 of file physics_model.h.