[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
periodic_1D_unsteady.cpp
1 #include "periodic_1D_unsteady.h"
2 
3 namespace PHiLiP {
4 
5 namespace FlowSolver {
6 
7 //=========================================================
8 // PERIODIC 1D DOMAIN FOR UNSTEADY CALCULATIONS
9 //=========================================================
10 
11 template <int dim, int nstate>
13  : PeriodicCubeFlow<dim, nstate>(parameters_input)
14  , unsteady_data_table_filename_with_extension(this->all_param.flow_solver_param.unsteady_data_table_filename+".txt")
15 {
16 
17 }
18 
19 
20 template <int dim, int nstate>
22  const std::shared_ptr <DGBase<dim, double>> dg
23  ) const
24 {
25  //Calculating energy via matrix-vector product
26  dealii::LinearAlgebra::distributed::Vector<double> temp;
27  temp.reinit(dg->solution);
29  dg->apply_global_mass_matrix(dg->solution, temp);
30  } else{
31  dg->global_mass_matrix.vmult(temp,dg->solution);
32  } //replace stage_j with M*stage_j
33  return temp * dg->solution;
34 }
35 
36 template <int dim, int nstate>
38  const std::shared_ptr <DGBase<dim, double>> dg
39  ) const
40 {
41  return compute_energy(dg);
42 }
43 
44 template <int dim, int nstate>
46  const unsigned int current_iteration,
47  const double current_time,
48  const std::shared_ptr <DGBase<dim, double>> dg ,
49  const std::shared_ptr <dealii::TableHandler> unsteady_data_table )
50 {
51  const double dt = this->all_param.ode_solver_param.initial_time_step;
52  int output_solution_every_n_iterations = round(this->all_param.ode_solver_param.output_solution_every_dt_time_intervals/dt);
53  if (this->all_param.ode_solver_param.output_solution_every_x_steps > output_solution_every_n_iterations)
54  output_solution_every_n_iterations = this->all_param.ode_solver_param.output_solution_every_x_steps;
55 
57  const PDEEnum pde_type = this->all_param.pde_type;
58 
59  if (pde_type == PDEEnum::advection){
60  if ((current_iteration % output_solution_every_n_iterations) == 0){
61  this->pcout << " Iter: " << current_iteration
62  << " Time: " << current_time
63  << std::endl;
64  }
65  (void) dg;
66  (void) unsteady_data_table;
67  }
68  else if (pde_type == PDEEnum::burgers_inviscid){
69  const double energy = this->compute_energy(dg);
70 
71  if ((current_iteration % output_solution_every_n_iterations) == 0){
72  this->pcout << " Iter: " << current_iteration
73  << " Time: " << current_time
74  << " Energy: " << energy
75  << std::endl;
76  }
77 
78  //detecting if the current run is calculating a reference solution
80  const double final_time = this->all_param.flow_solver_param.final_time;
81  const bool is_reference_solution = (dt < 2 * final_time/number_timesteps_ref);
82 
83  if(this->mpi_rank==0 && !is_reference_solution) {
84  //omit writing if current calculation is for a reference solution
85  unsteady_data_table->add_value("iteration", current_iteration);
86  this->add_value_to_data_table(current_time,"time",unsteady_data_table);
87  this->add_value_to_data_table(energy,"energy",unsteady_data_table);
88  std::ofstream unsteady_data_table_file(this->unsteady_data_table_filename_with_extension);
89  unsteady_data_table->write_text(unsteady_data_table_file);
90  }
91  }
92 
93 }
94 
95 #if PHILIP_DIM==1
97 #endif
98 
99 } // FlowSolver namespace
100 } // PHiLiP namespace
101 
const Parameters::AllParameters all_param
All parameters.
double final_time
Final solution time.
PartialDifferentialEquation pde_type
Store the PDE type to be solved.
double compute_energy(const std::shared_ptr< DGBase< dim, double >> dg) const
Calculate energy as a matrix-vector product, solution^T (M+K) solution.
TimeRefinementStudyParam time_refinement_study_param
Contains the parameters for time refinement study.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
Periodic1DUnsteady(const Parameters::AllParameters *const parameters_input)
Constructor.
int output_solution_every_x_steps
Outputs the solution every x steps to .vtk file.
PartialDifferentialEquation
Possible Partial Differential Equations to solve.
Files for the baseline physics.
Definition: ADTypes.hpp:10
int number_of_timesteps_for_reference_solution
For time refinement study with reference solution, number of steps for reference solution.
double output_solution_every_dt_time_intervals
Outputs the solution every dt time intervals to .vtk file.
Main parameter class that contains the various other sub-parameter classes.
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.
double initial_time_step
Time step used in ODE solver.
std::string unsteady_data_table_filename_with_extension
Filename for unsteady data.
bool use_inverse_mass_on_the_fly
Flag to use inverse mass matrix on-the-fly for explicit solves.
double get_numerical_entropy(const std::shared_ptr< DGBase< dim, double >> dg) const
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
dealii::ConditionalOStream pcout
ConditionalOStream.
void compute_unsteady_data_and_write_to_table(const unsigned int current_iteration, const double current_time, 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.