[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
homogeneous_isotropic_turbulence_initialization_check.cpp
1 #include "homogeneous_isotropic_turbulence_initialization_check.h"
2 #include "flow_solver/flow_solver_factory.h"
3 #include "flow_solver/flow_solver_cases/periodic_turbulence.h"
4 #include "physics/initial_conditions/set_initial_condition.h"
5 #include <deal.II/base/table_handler.h>
6 #include <algorithm>
7 #include <iterator>
8 #include <string>
9 #include <fstream>
10 
11 namespace PHiLiP {
12 namespace Tests {
13 
14 template <int dim, int nstate>
16  const PHiLiP::Parameters::AllParameters *const parameters_input,
17  const dealii::ParameterHandler &parameter_handler_input)
18  : TestsBase::TestsBase(parameters_input)
19  , parameter_handler(parameter_handler_input)
20  , kinetic_energy_expected(parameters_input->flow_solver_param.expected_kinetic_energy_at_final_time)
21 {}
22 
23 template <int dim, int nstate>
25 {
26  // copy all parameters
28 
29  std::unique_ptr<FlowSolver::FlowSolver<dim,nstate>> flow_solver = FlowSolver::FlowSolverFactory<dim,nstate>::select_flow_case(&parameters, parameter_handler);
30 
31  this->pcout << "Outputting solution files at initialization... " << std::flush;
32  flow_solver->dg->output_results_vtk(9999);
33  this->pcout << "done." << std::endl;
34 
35  std::unique_ptr<FlowSolver::PeriodicTurbulence<dim, nstate>> flow_solver_case = std::make_unique<FlowSolver::PeriodicTurbulence<dim,nstate>>(this->all_parameters);
36  flow_solver_case->output_velocity_field(flow_solver->dg,0,0.0);
37 
38  return 0;
39 }
40 
41 #if PHILIP_DIM==3
43 #endif
44 } // Tests namespace
45 } // PHiLiP namespace
Files for the baseline physics.
Definition: ADTypes.hpp:10
HomogeneousIsotropicTurbulenceInitializationCheck(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.
Main parameter class that contains the various other sub-parameter classes.
const dealii::ParameterHandler & parameter_handler
Parameter handler for storing the .prm file being ran.
static std::unique_ptr< FlowSolver< dim, nstate > > select_flow_case(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Factory to return the correct flow solver given input file.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
Definition: tests.h:20
dealii::ConditionalOStream pcout
ConditionalOStream.
Definition: tests.h:45
Base class of all the tests.
Definition: tests.h:17