1 #include "periodic_turbulence.h" 2 #include "ode_solver/low_storage_runge_kutta_ode_solver.h" 3 #include <deal.II/base/function.h> 6 #include <deal.II/dofs/dof_tools.h> 7 #include <deal.II/grid/grid_tools.h> 8 #include <deal.II/numerics/vector_tools.h> 9 #include <deal.II/fe/fe_values.h> 10 #include "physics/physics_factory.h" 11 #include <deal.II/base/table_handler.h> 12 #include <deal.II/base/tensor.h> 15 #include <deal.II/base/quadrature_lib.h> 19 namespace FlowSolver {
24 template <
int dim,
int nstate>
27 , unsteady_data_table_filename_with_extension(this->all_param.flow_solver_param.unsteady_data_table_filename+
".txt")
28 , number_of_times_to_output_velocity_field(this->all_param.flow_solver_param.number_of_times_to_output_velocity_field)
29 , output_velocity_field_at_fixed_times(this->all_param.flow_solver_param.output_velocity_field_at_fixed_times)
30 , output_vorticity_magnitude_field_in_addition_to_velocity(this->all_param.flow_solver_param.output_vorticity_magnitude_field_in_addition_to_velocity)
31 , output_flow_field_files_directory_name(this->all_param.flow_solver_param.output_flow_field_files_directory_name)
32 , output_solution_at_exact_fixed_times(this->all_param.ode_solver_param.output_solution_at_exact_fixed_times)
46 parameters_navier_stokes.
pde_type = Parameters::AllParameters::PartialDifferentialEquation::navier_stokes;
63 std::string line = output_velocity_field_times_string;
64 std::string::size_type sz1;
67 line = line.substr(sz1);
91 template <
int dim,
int nstate>
95 std::string flow_type_string;
103 template <
int dim,
int nstate>
108 this->
pcout <<
"- - Using constant time step in FlowSolver parameters: " << constant_time_step << std::endl;
109 return constant_time_step;
111 const unsigned int number_of_degrees_of_freedom_per_state = dg->dof_handler.n_dofs()/nstate;
112 const double approximate_grid_spacing = (this->
domain_right-this->
domain_left)/pow(number_of_degrees_of_freedom_per_state,(1.0/dim));
114 return constant_time_step;
118 std::string get_padded_mpi_rank_string(
const int mpi_rank_input) {
120 std::string mpi_rank_string = std::to_string(mpi_rank_input);
121 const unsigned int length_of_mpi_rank_with_padding = 5;
122 const int number_of_zeros = length_of_mpi_rank_with_padding - mpi_rank_string.length();
123 mpi_rank_string.insert(0, number_of_zeros,
'0');
125 return mpi_rank_string;
128 template<
int dim,
int nstate>
131 const unsigned int output_file_index,
132 const double current_time)
const 134 this->
pcout <<
" ... Writting velocity field ... " << std::flush;
144 const std::string mpi_rank_string = get_padded_mpi_rank_string(this->
mpi_rank);
146 const std::string filename_without_extension = filename_prefix + std::string(
"-") + mpi_rank_string;
158 std::ofstream data_table_file(filename_for_time_table);
165 std::ofstream FILE (filename);
168 if (!FILE.is_open()) {
169 this->
pcout <<
"ERROR: Cannot open file " << filename << std::endl;
172 const unsigned int number_of_degrees_of_freedom_per_state = dg->dof_handler.n_dofs()/nstate;
173 FILE << number_of_degrees_of_freedom_per_state << std::string(
"\n");
177 dealii::Quadrature<1> vol_quad_equidistant_1D = dealii::QIterated<1>(dealii::QTrapez<1>(),dg->max_degree);
178 dealii::FE_DGQArbitraryNodes<1,1> equidistant_finite_element(vol_quad_equidistant_1D);
180 const unsigned int init_grid_degree = dg->high_order_grid->fe_system.tensor_degree();
190 const unsigned int max_dofs_per_cell = dg->dof_handler.get_fe_collection().max_dofs_per_cell();
191 std::vector<dealii::types::global_dof_index> current_dofs_indices(max_dofs_per_cell);
192 auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
193 for (
auto current_cell = dg->dof_handler.begin_active(); current_cell!=dg->dof_handler.end(); ++current_cell, ++metric_cell) {
194 if (!current_cell->is_locally_owned())
continue;
196 const int i_fele = current_cell->active_fe_index();
197 const unsigned int poly_degree = i_fele;
198 const unsigned int n_dofs_cell = dg->fe_collection[poly_degree].dofs_per_cell;
199 const unsigned int n_shape_fns = n_dofs_cell / nstate;
200 const unsigned int n_quad_pts = n_shape_fns;
203 const dealii::FESystem<dim> &fe_metric = dg->high_order_grid->fe_system;
204 const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
205 const unsigned int n_grid_nodes = n_metric_dofs / dim;
206 std::vector<dealii::types::global_dof_index> metric_dof_indices(n_metric_dofs);
207 metric_cell->get_dof_indices (metric_dof_indices);
208 std::array<std::vector<double>,dim> mapping_support_points;
209 for(
int idim=0; idim<dim; idim++){
210 mapping_support_points[idim].resize(n_grid_nodes);
214 const std::vector<unsigned int > &index_renumbering = dealii::FETools::hierarchic_to_lexicographic_numbering<dim>(init_grid_degree);
215 for (
unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
216 const double val = (dg->high_order_grid->volume_nodes[metric_dof_indices[idof]]);
217 const unsigned int istate = fe_metric.system_to_component_index(idof).first;
218 const unsigned int ishape = fe_metric.system_to_component_index(idof).second;
219 const unsigned int igrid_node = index_renumbering[ishape];
220 mapping_support_points[istate][igrid_node] = val;
228 n_quad_pts, n_grid_nodes,
229 mapping_support_points,
230 mapping_basis_at_equidistant,
231 dg->all_parameters->use_invariant_curl_form);
234 current_dofs_indices.resize(n_dofs_cell);
235 current_cell->get_dof_indices (current_dofs_indices);
237 std::array<std::vector<double>,nstate> soln_coeff;
238 for(
unsigned int idof=0; idof<n_dofs_cell; idof++){
239 const unsigned int istate = dg->fe_collection[poly_degree].system_to_component_index(idof).first;
240 const unsigned int ishape = dg->fe_collection[poly_degree].system_to_component_index(idof).second;
242 soln_coeff[istate].resize(n_shape_fns);
244 soln_coeff[istate][ishape] = dg->solution(current_dofs_indices[idof]);
247 std::array<std::vector<double>,nstate> soln_at_q;
248 std::array<dealii::Tensor<1,dim,std::vector<double>>,nstate> soln_grad_at_q;
249 for(
int istate=0; istate<nstate; istate++){
250 soln_at_q[istate].resize(n_quad_pts);
255 dealii::Tensor<1,dim,std::vector<double>> ref_gradient_basis_fns_times_soln;
256 for(
int idim=0; idim<dim; idim++){
257 ref_gradient_basis_fns_times_soln[idim].resize(n_quad_pts);
264 for(
int idim=0; idim<dim; idim++){
265 soln_grad_at_q[istate][idim].resize(n_quad_pts);
266 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
267 for(
int jdim=0; jdim<dim; jdim++){
269 soln_grad_at_q[istate][idim][iquad] += metric_oper_equid.
metric_cofactor_vol[idim][jdim][iquad]
270 * ref_gradient_basis_fns_times_soln[jdim][iquad]
277 dealii::Tensor<1,dim,std::vector<double>> velocity_at_q;
278 std::vector<double> vorticity_magnitude_at_q(n_quad_pts);
279 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
280 std::array<double,nstate> soln_state;
281 std::array<dealii::Tensor<1,dim,double>,nstate> soln_grad_state;
282 for(
int istate=0; istate<nstate; istate++){
283 soln_state[istate] = soln_at_q[istate][iquad];
284 for(
int idim=0; idim<dim; idim++){
285 soln_grad_state[istate][idim] = soln_grad_at_q[istate][idim][iquad];
288 const dealii::Tensor<1,dim,double> velocity = this->
navier_stokes_physics->compute_velocities(soln_state);
289 for(
int idim=0; idim<dim; idim++){
291 velocity_at_q[idim].resize(n_quad_pts);
292 velocity_at_q[idim][iquad] = velocity[idim];
297 vorticity_magnitude_at_q[iquad] = this->
navier_stokes_physics->compute_vorticity_magnitude(soln_state, soln_grad_state);
301 for(
unsigned int ishape=0; ishape<n_shape_fns; ishape++){
302 dealii::Point<dim,double> vol_equid_node;
304 for(
int idim=0; idim<dim; idim++) {
305 vol_equid_node[idim] = metric_oper_equid.
flux_nodes_vol[idim][ishape];
306 FILE << std::setprecision(17) << vol_equid_node[idim] << std::string(
" ");
309 for (
int d=0; d<dim; ++d) {
310 FILE << std::setprecision(17) << velocity_at_q[d][ishape] << std::string(
" ");
314 FILE << std::setprecision(17) << vorticity_magnitude_at_q[ishape] << std::string(
" ");
316 FILE << std::string(
"\n");
320 this->
pcout <<
"done." << std::endl;
323 template<
int dim,
int nstate>
326 std::array<double,NUMBER_OF_INTEGRATED_QUANTITIES> integral_values;
327 std::fill(integral_values.begin(), integral_values.end(), 0.0);
333 int overintegrate = 10;
336 dealii::QGauss<dim> quad_extra(dg.
max_degree+1+overintegrate);
337 dealii::QGauss<1> quad_extra_1D(dg.
max_degree+1+overintegrate);
339 const unsigned int n_quad_pts = quad_extra.size();
340 const unsigned int grid_degree = dg.
high_order_grid->fe_system.tensor_degree();
341 const unsigned int poly_degree = dg.
max_degree;
351 const std::vector<double> &quad_weights = quad_extra.get_weights();
354 const bool store_vol_flux_nodes =
false;
355 const bool store_surf_flux_nodes =
false;
357 const unsigned int n_dofs = dg.
fe_collection[poly_degree].n_dofs_per_cell();
358 const unsigned int n_shape_fns = n_dofs / nstate;
359 std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs);
360 auto metric_cell = dg.
high_order_grid->dof_handler_grid.begin_active();
363 if (!cell->is_locally_owned())
continue;
364 cell->get_dof_indices (dofs_indices);
367 const dealii::FESystem<dim> &fe_metric = dg.
high_order_grid->fe_system;
368 const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
369 const unsigned int n_grid_nodes = n_metric_dofs / dim;
370 std::vector<dealii::types::global_dof_index> metric_dof_indices(n_metric_dofs);
371 metric_cell->get_dof_indices (metric_dof_indices);
372 std::array<std::vector<double>,dim> mapping_support_points;
373 for(
int idim=0; idim<dim; idim++){
374 mapping_support_points[idim].resize(n_grid_nodes);
378 const std::vector<unsigned int > &index_renumbering = dealii::FETools::hierarchic_to_lexicographic_numbering<dim>(grid_degree);
379 for (
unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
380 const double val = (dg.
high_order_grid->volume_nodes[metric_dof_indices[idof]]);
381 const unsigned int istate = fe_metric.system_to_component_index(idof).first;
382 const unsigned int ishape = fe_metric.system_to_component_index(idof).second;
383 const unsigned int igrid_node = index_renumbering[ishape];
384 mapping_support_points[istate][igrid_node] = val;
392 n_quad_pts, n_grid_nodes,
393 mapping_support_points,
402 std::array<std::vector<double>,nstate> soln_coeff;
403 for (
unsigned int idof = 0; idof < n_dofs; ++idof) {
404 const unsigned int istate = dg.
fe_collection[poly_degree].system_to_component_index(idof).first;
405 const unsigned int ishape = dg.
fe_collection[poly_degree].system_to_component_index(idof).second;
407 soln_coeff[istate].resize(n_shape_fns);
410 soln_coeff[istate][ishape] = dg.
solution(dofs_indices[idof]);
414 std::array<std::vector<double>,nstate> soln_at_q_vect;
415 std::array<dealii::Tensor<1,dim,std::vector<double>>,nstate> soln_grad_at_q_vect;
416 for(
int istate=0; istate<nstate; istate++){
417 soln_at_q_vect[istate].resize(n_quad_pts);
422 dealii::Tensor<1,dim,std::vector<double>> ref_gradient_basis_fns_times_soln;
423 for(
int idim=0; idim<dim; idim++){
424 ref_gradient_basis_fns_times_soln[idim].resize(n_quad_pts);
425 soln_grad_at_q_vect[istate][idim].resize(n_quad_pts);
432 for(
int idim=0; idim<dim; idim++){
433 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
434 for(
int jdim=0; jdim<dim; jdim++){
436 soln_grad_at_q_vect[istate][idim][iquad] += metric_oper.
metric_cofactor_vol[idim][jdim][iquad]
437 * ref_gradient_basis_fns_times_soln[jdim][iquad]
445 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
447 std::array<double,nstate> soln_at_q;
448 std::array<dealii::Tensor<1,dim,double>,nstate> soln_grad_at_q;
450 for(
int istate=0; istate<nstate; istate++){
451 soln_at_q[istate] = soln_at_q_vect[istate][iquad];
452 for(
int idim=0; idim<dim; idim++){
453 soln_grad_at_q[istate][idim] = soln_grad_at_q_vect[istate][idim][iquad];
457 std::array<double,NUMBER_OF_INTEGRATED_QUANTITIES> integrand_values;
458 std::fill(integrand_values.begin(), integrand_values.end(), 0.0);
459 integrand_values[IntegratedQuantitiesEnum::kinetic_energy] = this->
navier_stokes_physics->compute_kinetic_energy_from_conservative_solution(soln_at_q);
460 integrand_values[IntegratedQuantitiesEnum::enstrophy] = this->
navier_stokes_physics->compute_enstrophy(soln_at_q,soln_grad_at_q);
461 integrand_values[IntegratedQuantitiesEnum::pressure_dilatation] = this->
navier_stokes_physics->compute_pressure_dilatation(soln_at_q,soln_grad_at_q);
462 integrand_values[IntegratedQuantitiesEnum::deviatoric_strain_rate_tensor_magnitude_sqr] = this->
navier_stokes_physics->compute_deviatoric_strain_rate_tensor_magnitude_sqr(soln_at_q,soln_grad_at_q);
463 integrand_values[IntegratedQuantitiesEnum::strain_rate_tensor_magnitude_sqr] = this->
navier_stokes_physics->compute_strain_rate_tensor_magnitude_sqr(soln_at_q,soln_grad_at_q);
466 integral_values[i_quantity] += integrand_values[i_quantity] * quad_weights[iquad] * metric_oper.
det_Jac_vol[iquad];
485 template<
int dim,
int nstate>
488 const double integrated_kinetic_energy = this->
integrated_quantities[IntegratedQuantitiesEnum::kinetic_energy];
495 return integrated_kinetic_energy;
498 template<
int dim,
int nstate>
504 template<
int dim,
int nstate>
507 const double integrated_enstrophy = this->
integrated_quantities[IntegratedQuantitiesEnum::enstrophy];
508 double vorticity_based_dissipation_rate = 0.0;
510 vorticity_based_dissipation_rate = this->
navier_stokes_physics->compute_vorticity_based_dissipation_rate_from_integrated_enstrophy(integrated_enstrophy);
512 return vorticity_based_dissipation_rate;
515 template<
int dim,
int nstate>
518 const double integrated_pressure_dilatation = this->
integrated_quantities[IntegratedQuantitiesEnum::pressure_dilatation];
519 return (-1.0*integrated_pressure_dilatation);
522 template<
int dim,
int nstate>
525 const double integrated_deviatoric_strain_rate_tensor_magnitude_sqr = this->
integrated_quantities[IntegratedQuantitiesEnum::deviatoric_strain_rate_tensor_magnitude_sqr];
526 double deviatoric_strain_rate_tensor_based_dissipation_rate = 0.0;
528 deviatoric_strain_rate_tensor_based_dissipation_rate =
529 this->
navier_stokes_physics->compute_deviatoric_strain_rate_tensor_based_dissipation_rate_from_integrated_deviatoric_strain_rate_tensor_magnitude_sqr(integrated_deviatoric_strain_rate_tensor_magnitude_sqr);
531 return deviatoric_strain_rate_tensor_based_dissipation_rate;
534 template<
int dim,
int nstate>
537 const double integrated_strain_rate_tensor_magnitude_sqr = this->
integrated_quantities[IntegratedQuantitiesEnum::strain_rate_tensor_magnitude_sqr];
538 double strain_rate_tensor_based_dissipation_rate = 0.0;
540 strain_rate_tensor_based_dissipation_rate =
541 this->
navier_stokes_physics->compute_strain_rate_tensor_based_dissipation_rate_from_integrated_strain_rate_tensor_magnitude_sqr(integrated_strain_rate_tensor_magnitude_sqr);
543 return strain_rate_tensor_based_dissipation_rate;
546 template<
int dim,
int nstate>
552 template<
int dim,
int nstate>
559 const unsigned int n_dofs_cell = dg->fe_collection[poly_degree].dofs_per_cell;
560 const unsigned int n_quad_pts = dg->volume_quadrature_collection[poly_degree].size();
561 const unsigned int n_shape_fns = n_dofs_cell / nstate;
564 vol_projection.
build_1D_volume_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
568 soln_basis.
build_1D_volume_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
574 std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs_cell);
576 double integrand_numerical_entropy_function=0;
577 double integral_numerical_entropy_function=0;
578 const std::vector<double> &quad_weights = dg->volume_quadrature_collection[poly_degree].get_weights();
580 auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
582 for (
auto cell = dg->dof_handler.begin_active(); cell!= dg->dof_handler.end(); ++cell, ++metric_cell) {
583 if (!cell->is_locally_owned())
continue;
584 cell->get_dof_indices (dofs_indices);
587 const dealii::FESystem<dim> &fe_metric = dg->high_order_grid->fe_system;
588 const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
589 const unsigned int n_grid_nodes = n_metric_dofs / dim;
590 std::vector<dealii::types::global_dof_index> metric_dof_indices(n_metric_dofs);
591 metric_cell->get_dof_indices (metric_dof_indices);
592 std::array<std::vector<double>,dim> mapping_support_points;
593 for(
int idim=0; idim<dim; idim++){
594 mapping_support_points[idim].resize(n_grid_nodes);
598 const std::vector<unsigned int > &index_renumbering = dealii::FETools::hierarchic_to_lexicographic_numbering<dim>(dg->max_grid_degree);
599 for (
unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
600 const double val = (dg->high_order_grid->volume_nodes[metric_dof_indices[idof]]);
601 const unsigned int istate = fe_metric.system_to_component_index(idof).first;
602 const unsigned int ishape = fe_metric.system_to_component_index(idof).second;
603 const unsigned int igrid_node = index_renumbering[ishape];
604 mapping_support_points[istate][igrid_node] = val;
612 n_quad_pts, n_grid_nodes,
613 mapping_support_points,
615 dg->all_parameters->use_invariant_curl_form);
622 std::array<std::vector<double>,nstate> soln_coeff;
623 for (
unsigned int idof = 0; idof < n_dofs_cell; ++idof) {
624 const unsigned int istate = dg->fe_collection[poly_degree].system_to_component_index(idof).first;
625 const unsigned int ishape = dg->fe_collection[poly_degree].system_to_component_index(idof).second;
627 soln_coeff[istate].resize(n_shape_fns);
629 soln_coeff[istate][ishape] = dg->solution(dofs_indices[idof]);
633 std::array<std::vector<double>,nstate> soln_at_q;
634 for(
int istate=0; istate<nstate; istate++){
635 soln_at_q[istate].resize(n_quad_pts);
642 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
644 std::array<double,nstate> soln_state;
646 for(
int istate=0; istate<nstate; istate++){
647 soln_state[istate] = soln_at_q[istate][iquad];
649 integrand_numerical_entropy_function = this->
navier_stokes_physics->compute_numerical_entropy_function(soln_state);
650 integral_numerical_entropy_function += integrand_numerical_entropy_function * quad_weights[iquad] * metric_oper.
det_Jac_vol[iquad];
654 const double mpi_integrated_numerical_entropy = dealii::Utilities::MPI::sum(integral_numerical_entropy_function, this->
mpi_communicator);
656 return mpi_integrated_numerical_entropy;
660 template <
int dim,
int nstate>
662 const double FR_entropy_contribution_RRK_solver,
663 const unsigned int current_iteration,
669 if (current_iteration==0) {
680 template <
int dim,
int nstate>
684 const std::shared_ptr <dealii::TableHandler> unsteady_data_table)
687 const unsigned int current_iteration = ode_solver->current_iteration;
688 const double current_time = ode_solver->current_time;
701 const double relaxation_parameter = ode_solver->relaxation_parameter_RRK_solver;
720 unsteady_data_table->write_text(unsteady_data_table_file);
723 this->
pcout <<
" Iter: " << current_iteration
724 <<
" Time: " << current_time
725 <<
" Energy: " << integrated_kinetic_energy
726 <<
" Enstrophy: " << integrated_enstrophy;
728 this->
pcout <<
" eps_vorticity: " << vorticity_based_dissipation_rate
729 <<
" eps_p+eps_strain: " << (pressure_dilatation_based_dissipation_rate + strain_rate_tensor_based_dissipation_rate);
735 this->
pcout <<
" Relaxation Parameter: " << std::setprecision(16) << relaxation_parameter;
737 this->
pcout << std::endl;
740 if(std::isnan(integrated_kinetic_energy)) {
741 this->
pcout <<
" ERROR: Kinetic energy at time " << current_time <<
" is nan." << std::endl;
742 this->
pcout <<
" Consider decreasing the time step / CFL number." << std::endl;
749 const double next_time = current_time +
time_step;
752 bool is_output_time =
false;
754 is_output_time = current_time == desired_time;
756 is_output_time = ((current_time<=desired_time) && (next_time>desired_time));
const Parameters::AllParameters all_param
All parameters.
const MPI_Comm mpi_communicator
MPI communicator.
double get_integrated_enstrophy() const
FlowCaseType
Selects the flow case to be simulated.
dealii::Tensor< 2, dim, std::vector< real > > metric_cofactor_vol
The volume metric cofactor matrix.
PartialDifferentialEquation pde_type
Store the PDE type to be solved.
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
double initial_time
Initial time at which we initialize the ODE solver with.
double get_vorticity_based_dissipation_rate() const
FlowCaseType flow_case_type
Selected FlowCaseType from the input file.
double courant_friedrichs_lewy_number
Courant-Friedrich-Lewy (CFL) number for constant time step.
bool adaptive_time_step
Flag for computing the time step on the fly.
const bool output_vorticity_magnitude_field_in_addition_to_velocity
Flag for outputting vorticity magnitude field in addition to velocity field at fixed times...
const bool output_solution_at_exact_fixed_times
Flag for outputting the solution at exact fixed times by decreasing the time step on the fly...
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
double constant_time_step
Constant time step.
dealii::FullMatrix< double > oneD_grad_operator
Stores the one dimensional gradient operator.
double mach_inf
Mach number at infinity.
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
double get_integrated_kinetic_energy() const
bool restart_computation_from_file
Restart computation from restart file.
double get_time_step() const
Getter for time step.
double get_strain_rate_tensor_based_dissipation_rate() const
double get_deviatoric_strain_rate_tensor_based_dissipation_rate() const
dealii::QGauss< 0 > oneD_face_quadrature
1D surface quadrature is always one single point for all poly degrees.
unsigned int index_of_current_desired_time_to_output_velocity_field
Index of current desired time to output velocity field.
Files for the baseline physics.
void output_velocity_field(std::shared_ptr< DGBase< dim, double >> dg, const unsigned int output_file_index, const double current_time) const
Output the velocity field to file.
const std::string output_flow_field_files_directory_name
Directory for writting flow field files.
const std::string unsteady_data_table_filename_with_extension
Filename (with extension) for the unsteady data table.
const double domain_left
Domain left-boundary value for generating the grid.
double initial_numerical_entropy_abs
Numerical entropy at initial time.
double reynolds_number_inf
Farfield Reynolds number.
double get_numerical_entropy(const std::shared_ptr< DGBase< dim, double >>) const
Retrieves cumulative_numerical_entropy_change_FRcorrected.
std::shared_ptr< dealii::TableHandler > exact_output_times_of_velocity_field_files_table
Data table storing the exact output times for the velocity field files.
EulerParam euler_param
Contains parameters for the Euler equations non-dimensionalization.
ODESolverEnum
Types of ODE solver.
double previous_numerical_entropy
Numerical entropy at previous timestep.
void build_1D_shape_functions_at_grid_nodes(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Constructs the volume operator and gradient operator.
unsigned int poly_degree
Polynomial order (P) of the basis functions for DG.
Main parameter class that contains the various other sub-parameter classes.
const dealii::hp::FECollection< 1 > oneD_fe_collection_1state
1D Finite Element Collection for p-finite-element to represent the solution for a single state...
bool do_calculate_numerical_entropy
For TGV, flag to calculate and write numerical entropy.
const bool output_velocity_field_at_fixed_times
Flag for outputting velocity field at fixed times.
bool do_calculate_numerical_entropy
Identifies if numerical entropy should be calculated; initialized as false.
std::string output_velocity_field_times_string
String of velocity field output times.
dealii::LinearAlgebra::distributed::Vector< double > solution
Current modal coefficients of the solution.
std::vector< real > det_Jac_vol
The determinant of the metric Jacobian at volume cubature nodes.
void display_additional_flow_case_specific_parameters() const override
Display additional more specific flow case parameters.
void add_value_to_data_table(const double value, const std::string value_string, const std::shared_ptr< dealii::TableHandler > data_table) const
Add a value to a given data table with scientific format.
ODESolverParam ode_solver_param
Contains parameters for ODE solver.
PeriodicTurbulence(const Parameters::AllParameters *const parameters_input)
Constructor.
dealii::Table< 1, double > output_velocity_field_times
Times at which to output the velocity field.
NavierStokesParam navier_stokes_param
Contains parameters for the Navier-Stokes equations non-dimensionalization.
void build_1D_gradient_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
bool is_viscous_flow
Identifies if viscous flow; initialized as true.
const unsigned int max_degree
Maximum degree used for p-refi1nement.
Base metric operators class that stores functions used in both the volume and on surface.
void compute_and_update_integrated_quantities(DGBase< dim, double > &dg)
double get_constant_time_step(std::shared_ptr< DGBase< dim, double >> dg) const override
Function to compute the constant time step.
The mapping shape functions evaluated at the desired nodes (facet set included in volume grid nodes f...
const unsigned int number_of_times_to_output_velocity_field
Number of times to output the velocity field.
bool is_decaying_homogeneous_isotropic_turbulence
Identified if DHIT case; initialized as false.
std::shared_ptr< HighOrderGrid< dim, real, MeshType > > high_order_grid
High order grid that will provide the MappingFEField.
void gradient_matrix_vector_mult_1D(const std::vector< real > &input_vect, dealii::Tensor< 1, dim, std::vector< real >> &output_vect, const dealii::FullMatrix< double > &basis, const dealii::FullMatrix< double > &gradient_basis)
Computes the gradient of a scalar using sum-factorization where the basis are the same in each direct...
void compute_unsteady_data_and_write_to_table(const std::shared_ptr< ODE::ODESolverBase< dim, double >> ode_solver, const std::shared_ptr< DGBase< dim, double >> dg, const std::shared_ptr< dealii::TableHandler > unsteady_data_table) override
Compute the desired unsteady data and write it to a table.
bool is_taylor_green_vortex
Identifies if taylor green vortex case; initialized as false.
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
double cumulative_numerical_entropy_change_FRcorrected
Cumulative change in numerical entropy.
static std::shared_ptr< PhysicsBase< dim, nstate, real > > create_Physics(const Parameters::AllParameters *const parameters_input, std::shared_ptr< ModelBase< dim, nstate, real > > model_input=nullptr)
Factory to return the correct physics given input file.
double compute_current_integrated_numerical_entropy(const std::shared_ptr< DGBase< dim, double >> dg) const
Calculate numerical entropy by matrix-vector product.
ODESolverEnum ode_solver_type
ODE solver type.
dealii::DoFHandler< dim > dof_handler
Finite Element Collection to represent the high-order grid.
std::shared_ptr< Physics::NavierStokes< dim, dim+2, double > > navier_stokes_physics
Pointer to Navier-Stokes physics object for computing things on the fly.
const double domain_size
Domain size (length in 1D, area in 2D, and volume in 3D)
void build_volume_metric_operators(const unsigned int n_quad_pts, const unsigned int n_metric_dofs, const std::array< std::vector< real >, dim > &mapping_support_points, mapping_shape_functions< dim, n_faces, real > &mapping_basis, const bool use_invariant_curl_form=false)
Builds the volume metric operators.
const double domain_right
Domain right-boundary value for generating the grid.
double maximum_local_wave_speed
Maximum local wave speed (i.e. convective eigenvalue)
std::array< double, NUMBER_OF_INTEGRATED_QUANTITIES > integrated_quantities
Array for storing the integrated quantities; done for computational efficiency.
const int mpi_rank
MPI rank.
std::string flow_field_quantity_filename_prefix
Flow field quantity filename prefix.
DGBase is independent of the number of state variables.
dealii::ConditionalOStream pcout
ConditionalOStream.
void build_1D_shape_functions_at_flux_nodes(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature, const dealii::Quadrature< 0 > &face_quadrature)
Constructs the volume, gradient, surface, and surface gradient operator.
double get_pressure_dilatation_based_dissipation_rate() const
void update_numerical_entropy(const double FR_entropy_contribution_RRK_solver, const unsigned int current_iteration, const std::shared_ptr< DGBase< dim, double >> dg)
Update numerical entropy variables.
void matrix_vector_mult_1D(const std::vector< real > &input_vect, std::vector< real > &output_vect, const dealii::FullMatrix< double > &basis_x, const bool adding=false, const double factor=1.0)
Apply the matrix vector operation using the 1D operator in each direction.
const dealii::hp::FECollection< dim > fe_collection
Finite Element Collection for p-finite-element to represent the solution.
dealii::Tensor< 1, dim, std::vector< real > > flux_nodes_vol
Stores the physical volume flux nodes.
Navier-Stokes equations. Derived from Euler for the convective terms, which is derived from PhysicsBa...
Projection operator corresponding to basis functions onto M-norm (L2).
void display_grid_parameters() const
Display grid parameters.
double time_step
Current time step.
static const int NUMBER_OF_INTEGRATED_QUANTITIES