4 #include <deal.II/base/tensor.h>     6 #include "parameters/all_parameters.h"     7 #include "parameters/parameters_manufactured_solution.h"    29 template <
int dim, 
int nstate, 
typename real>
    58         const bool                                                convection = 
true, 
    59         const bool                                                diffusion = 
true,
    66     std::array<dealii::Tensor<1,dim,real>,nstate> 
convective_flux (
const std::array<real,nstate> &solution) 
const;
    70                 const std::array<real,nstate> &conservative_soln1,
    71                 const std::array<real,nstate> &conservative_soln2) 
const override;
    75                 const std::array<real,nstate> &conservative_soln) 
const;
    79                 const std::array<real,nstate> &entropy_var) 
const;
    83         const std::array<real,nstate> &,
    84         const dealii::Tensor<1,dim,real> &) 
const;
    99         const std::array<real,nstate> &solution,
   100         const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
   101         const dealii::types::global_dof_index cell_index) 
const;
   105         const std::array<real,nstate> &solution,
   106         const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient) 
const;
   110         const dealii::Point<dim,real> &pos,
   111         const std::array<real,nstate> &solution,
   112         const real current_time,
   113         const dealii::types::global_dof_index cell_index) 
const;
   117         const dealii::Point<dim,real> &pos,
   118         const std::array<real,nstate> &solution,
   119         const real current_time) 
const;
   127         const dealii::Point<dim, real> &,
   128         const dealii::Tensor<1,dim,real> &,
   129         const std::array<real,nstate> &,
   130         const std::array<dealii::Tensor<1,dim,real>,nstate> &,
   131         std::array<real,nstate> &,
   132         std::array<dealii::Tensor<1,dim,real>,nstate> &) 
const;
 real max_viscous_eigenvalue(const std::array< real, nstate > &soln) const
Maximum viscous eigenvalue. 
std::array< real, nstate > compute_entropy_variables(const std::array< real, nstate > &conservative_soln) const
Computes the entropy variables. 
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 override
Convective split flux. 
TestType
Possible integration tests to run. 
std::array< dealii::Tensor< 1, dim, real >, nstate > convective_flux(const std::array< real, nstate > &solution) const
Convective flux: . 
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived. 
Manufactured solution used for grid studies to check convergence orders. 
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
Dissipative flux: u. 
Files for the baseline physics. 
Main parameter class that contains the various other sub-parameter classes. 
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
Spectral radius of convective term Jacobian is 'c'. 
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
If diffusion is present, assign Dirichlet boundary condition. 
real max_convective_eigenvalue(const std::array< real, nstate > &soln) const
Maximum convective eigenvalue. 
static dealii::Tensor< 2, 3, double > get_default_diffusion_tensor()
gets the default diffusion tensor 
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
Source term is zero or depends on manufactured solution. 
double diffusion_scaling_coeff
Diffusion scaling coefficient in front of the diffusion tensor. 
const bool has_nonzero_physical_source
Flag to signal that physical source term is non-zero. 
const Parameters::AllParameters::TestType test_type
Allows Burgers to distinguish between different unsteady test types. 
Burgers(const Parameters::AllParameters *const parameters_input, const double diffusion_coefficient, const bool convection=true, const bool diffusion=true, const dealii::Tensor< 2, 3, double > input_diffusion_tensor=Parameters::ManufacturedSolutionParam::get_default_diffusion_tensor(), std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function=nullptr, const Parameters::AllParameters::TestType parameters_test=Parameters::AllParameters::TestType::run_control, const bool has_nonzero_physical_source=false)
Constructor. 
Burger's equation with nonlinear advective term and linear diffusive term. Derived from PhysicsBase...
const bool hasDiffusion
Turns on diffusive part of the Burgers problem. 
std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function
Manufactured solution function. 
const bool hasConvection
Turns on convective part of the Burgers problem. 
real diffusion_coefficient() const
Diffusion coefficient.