1 #include "burgers_stability.h" 3 #include <deal.II/base/convergence_table.h> 4 #include <deal.II/base/function.h> 5 #include <deal.II/base/function_parser.h> 6 #include <deal.II/base/numbers.h> 7 #include <deal.II/base/tensor.h> 8 #include <deal.II/grid/grid_generator.h> 9 #include <deal.II/grid/grid_in.h> 10 #include <deal.II/grid/grid_out.h> 11 #include <deal.II/grid/grid_refinement.h> 12 #include <deal.II/grid/grid_tools.h> 13 #include <deal.II/numerics/data_out.h> 14 #include <deal.II/numerics/solution_transfer.h> 15 #include <deal.II/numerics/vector_tools.h> 19 #include "dg/dg_base.hpp" 20 #include "dg/dg_factory.hpp" 21 #include "ode_solver/ode_solver_base.h" 22 #include "ode_solver/ode_solver_factory.h" 23 #include "parameters/all_parameters.h" 24 #include "parameters/parameters.h" 25 #include "physics/initial_conditions/initial_condition_function.h" 26 #include "physics/initial_conditions/set_initial_condition.h" 31 template <
int dim,
int nstate>
36 template<
int dim,
int nstate>
40 dealii::LinearAlgebra::distributed::Vector<double> mass_matrix_times_solution(dg->right_hand_side);
41 if(dg->all_parameters->use_inverse_mass_on_the_fly)
42 dg->apply_global_mass_matrix(dg->solution,mass_matrix_times_solution);
44 dg->global_mass_matrix.vmult( mass_matrix_times_solution, dg->solution);
47 energy = dg->solution * mass_matrix_times_solution;
52 template<
int dim,
int nstate>
56 double conservation = 0.0;
57 dealii::LinearAlgebra::distributed::Vector<double> mass_matrix_times_solution(dg->right_hand_side);
58 if(dg->all_parameters->use_inverse_mass_on_the_fly)
59 dg->apply_global_mass_matrix(dg->solution,mass_matrix_times_solution);
61 dg->global_mass_matrix.vmult( mass_matrix_times_solution, dg->solution);
63 const unsigned int n_dofs_cell = dg->fe_collection[poly_degree].dofs_per_cell;
64 const unsigned int n_quad_pts = dg->volume_quadrature_collection[poly_degree].size();
65 std::vector<double> ones(n_quad_pts, 1.0);
67 std::vector<double> ones_hat(n_dofs_cell);
70 vol_projection.
build_1D_volume_operator(dg->oneD_fe_collection[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
74 dealii::LinearAlgebra::distributed::Vector<double> ones_hat_global(dg->right_hand_side);
75 std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs_cell);
76 for (
auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
77 if (!cell->is_locally_owned())
continue;
78 cell->get_dof_indices (dofs_indices);
79 for(
unsigned int idof=0;idof<n_dofs_cell; idof++){
80 ones_hat_global[dofs_indices[idof]] = ones_hat[idof];
84 conservation = ones_hat_global * mass_matrix_times_solution;
89 template <
int dim,
int nstate>
92 pcout <<
" Running Burgers energy stability. " << std::endl;
97 const unsigned int n_grids = (all_parameters_new.
use_energy) ? 4 : 5;
98 std::vector<double> grid_size(n_grids);
99 std::vector<double> soln_error(n_grids);
100 unsigned int poly_degree = 4;
101 dealii::ConvergenceTable convergence_table;
102 const unsigned int igrid_start = 3;
103 const unsigned int grid_degree = 1;
105 for(
unsigned int igrid = igrid_start; igrid<n_grids; igrid++){
107 #if PHILIP_DIM==1 // dealii::parallel::distributed::Triangulation<dim> does not work for 1D 108 using Triangulation = dealii::Triangulation<dim>;
109 std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
110 typename dealii::Triangulation<dim>::MeshSmoothing(
111 dealii::Triangulation<dim>::smoothing_on_refinement |
112 dealii::Triangulation<dim>::smoothing_on_coarsening));
114 using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
115 std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
117 typename dealii::Triangulation<dim>::MeshSmoothing(
118 dealii::Triangulation<dim>::smoothing_on_refinement |
119 dealii::Triangulation<dim>::smoothing_on_coarsening));
122 dealii::GridGenerator::hyper_cube(*grid, left, right,
true);
125 std::vector<dealii::GridTools::PeriodicFacePair<typename dealii::Triangulation<PHILIP_DIM>::cell_iterator> > matched_pairs;
126 dealii::GridTools::collect_periodic_faces(*grid,0,1,0,matched_pairs);
127 grid->add_periodicity(matched_pairs);
129 std::vector<dealii::GridTools::PeriodicFacePair<typename dealii::parallel::distributed::Triangulation<PHILIP_DIM>::cell_iterator> > matched_pairs;
130 dealii::GridTools::collect_periodic_faces(*grid,0,1,0,matched_pairs);
131 if(dim>=2) dealii::GridTools::collect_periodic_faces(*grid,2,3,1,matched_pairs);
132 if(dim>=3) dealii::GridTools::collect_periodic_faces(*grid,4,5,2,matched_pairs);
133 grid->add_periodicity(matched_pairs);
135 grid->refine_global(igrid);
136 pcout <<
"Grid generated and refined" << std::endl;
138 const unsigned int n_global_active_cells2 = grid->n_global_active_cells();
139 double n_dofs_cfl = pow(n_global_active_cells2,dim) * pow(poly_degree+1.0, dim);
140 double delta_x = (right-left)/pow(n_dofs_cfl,(1.0/dim));
147 pcout <<
"dg created" <<std::endl;
148 dg->allocate_system (
false,
false,
false);
151 pcout<<
"Setting up Initial Condition"<<std::endl;
153 std::shared_ptr< InitialConditionFunction<dim,nstate,double> > initial_condition_function =
160 double finalTime = 3.0;
167 ode_solver->current_iteration = 0;
168 ode_solver->advance_solution_time(0.000001);
176 std::ofstream myfile (
"energy_plot_burgers.gpl" , std::ios::trunc);
178 ode_solver->current_iteration = 0;
179 for (
int i = 0; i < std::ceil(finalTime/dt); ++ i)
181 ode_solver->advance_solution_time(dt);
184 current_energy /=initial_energy;
185 std::cout << std::setprecision(16) << std::fixed;
186 pcout <<
"Energy at time " << i * dt <<
" is " << current_energy << std::endl;
187 myfile << i * dt <<
" " << std::fixed << std::setprecision(16) << current_energy << std::endl;
188 if (current_energy*initial_energy - initial_energy >= 1.0)
190 pcout<<
"Energy not monotonicaly decreasing"<<std::endl;
194 if ( (current_energy*initial_energy - initial_energy >= 1.0e-11)&&(all_parameters_new.
conv_num_flux_type == Parameters::AllParameters::ConvectiveNumericalFlux::two_point_flux) )
196 pcout<<
"Energy not conserved"<<std::endl;
202 current_conservation /=initial_conservation;
203 std::cout << std::setprecision(16) << std::fixed;
204 pcout <<
"Normalized Conservation at time " << i * dt <<
" is " << current_conservation<< std::endl;
205 myfile << i * dt <<
" " << std::fixed << std::setprecision(16) << current_conservation << std::endl;
206 if (current_conservation*initial_conservation - initial_conservation >= 10.00)
208 pcout <<
"Not conserved" << std::endl;
216 std::ofstream myfile2 (
"solution_burgers.gpl" , std::ios::trunc);
218 dealii::QGaussLobatto<dim> quad_extra(dg->max_degree+1);
219 dealii::FEValues<dim,dim> fe_values_extra(*(dg->high_order_grid->mapping_fe_field), dg->fe_collection[poly_degree], quad_extra,
220 dealii::update_values | dealii::update_JxW_values | dealii::update_quadrature_points);
221 const unsigned int n_quad_pts = fe_values_extra.n_quadrature_points;
222 std::array<double,nstate> soln_at_q;
223 std::vector<dealii::types::global_dof_index> dofs_indices (fe_values_extra.dofs_per_cell);
225 for (
auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
226 if (!cell->is_locally_owned())
continue;
228 fe_values_extra.reinit (cell);
229 cell->get_dof_indices (dofs_indices);
231 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
232 std::fill(soln_at_q.begin(), soln_at_q.end(), 0.0);
233 for (
unsigned int idof=0; idof<fe_values_extra.dofs_per_cell; ++idof) {
234 const unsigned int istate = fe_values_extra.get_fe().system_to_component_index(idof).first;
235 soln_at_q[istate] += dg->solution[dofs_indices[idof]] * fe_values_extra.shape_value_component(idof, iquad, istate);
237 const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
239 std::cout << std::setprecision(16) << std::fixed;
240 myfile2<< std::fixed << std::setprecision(16) << qpoint[0] << std::fixed << std::setprecision(16) <<
" " << soln_at_q[0]<< std::endl;
248 ode_solver->current_iteration = 0;
250 ode_solver->advance_solution_time(finalTime);
251 const unsigned int n_global_active_cells = grid->n_global_active_cells();
252 const unsigned int n_dofs = dg->dof_handler.n_dofs();
253 pcout <<
"Dimension: " << dim
254 <<
"\t Polynomial degree p: " << poly_degree
256 <<
"Grid number: " << igrid+1 <<
"/" << n_grids
257 <<
". Number of active cells: " << n_global_active_cells
258 <<
". Number of degrees of freedom: " << n_dofs
262 int overintegrate = 10;
263 dealii::QGauss<dim> quad_extra(poly_degree+1+overintegrate);
264 dealii::FEValues<dim,dim> fe_values_extra(*(dg->high_order_grid->mapping_fe_field), dg->fe_collection[poly_degree], quad_extra,
265 dealii::update_values | dealii::update_JxW_values | dealii::update_quadrature_points);
266 const unsigned int n_quad_pts = fe_values_extra.n_quadrature_points;
267 std::array<double,nstate> soln_at_q;
269 double l2error = 0.0;
272 const double pi = atan(1)*4.0;
273 std::vector<dealii::types::global_dof_index> dofs_indices (fe_values_extra.dofs_per_cell);
274 for (
auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
276 if (!cell->is_locally_owned())
continue;
278 fe_values_extra.reinit (cell);
279 cell->get_dof_indices (dofs_indices);
281 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
283 std::fill(soln_at_q.begin(), soln_at_q.end(), 0.0);
284 for (
unsigned int idof=0; idof<fe_values_extra.dofs_per_cell; ++idof) {
285 const unsigned int istate = fe_values_extra.get_fe().system_to_component_index(idof).first;
286 soln_at_q[istate] += dg->solution[dofs_indices[idof]] * fe_values_extra.shape_value_component(idof, iquad, istate);
289 for (
int istate=0; istate<nstate; ++istate) {
290 const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
292 for(
int idim=0; idim<dim; idim++){
293 uexact += cos(pi*(qpoint[idim]-finalTime));
295 l2error += pow(soln_at_q[istate] - uexact, 2) * fe_values_extra.JxW(iquad);
299 const double l2error_mpi_sum = std::sqrt(dealii::Utilities::MPI::sum(l2error,
mpi_communicator));
302 const double dx = 1.0/pow(n_dofs,(1.0/dim));
303 grid_size[igrid] = dx;
304 soln_error[igrid] = l2error_mpi_sum;
306 convergence_table.add_value(
"p", poly_degree);
307 convergence_table.add_value(
"cells", n_global_active_cells);
308 convergence_table.add_value(
"DoFs", n_dofs);
309 convergence_table.add_value(
"dx", dx);
310 convergence_table.add_value(
"soln_L2_error", l2error_mpi_sum);
312 pcout <<
" Grid size h: " << dx
313 <<
" L2-soln_error: " << l2error_mpi_sum
314 <<
" Residual: " << ode_solver->residual_norm
317 if (igrid > igrid_start) {
318 const double slope_soln_err = log(soln_error[igrid]/soln_error[igrid-1])
319 / log(grid_size[igrid]/grid_size[igrid-1]);
320 pcout <<
"From grid " << igrid
321 <<
" to grid " << igrid+1
322 <<
" dimension: " << dim
323 <<
" polynomial degree p: " << poly_degree
325 <<
" solution_error1 " << soln_error[igrid-1]
326 <<
" solution_error2 " << soln_error[igrid]
327 <<
" slope " << slope_soln_err
329 if(igrid == n_grids-1){
330 if(std::abs(slope_soln_err-(poly_degree+1))>0.05){
336 pcout <<
" ********************************************" 338 <<
" Convergence rates for p = " << poly_degree
340 <<
" ********************************************" 342 convergence_table.evaluate_convergence_rates(
"soln_L2_error",
"cells", dealii::ConvergenceTable::reduction_rate_log2, dim);
343 convergence_table.set_scientific(
"dx",
true);
344 convergence_table.set_scientific(
"soln_L2_error",
true);
345 if (
pcout.is_active()) convergence_table.write_text(
pcout.get_stream());
double compute_energy(std::shared_ptr< PHiLiP::DGBase< dim, double > > &dg) const
Function computes the energy.
Burgers' periodic unsteady test.
const MPI_Comm mpi_communicator
MPI communicator.
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
bool use_energy
Flag to use an energy monotonicity test.
Files for the baseline physics.
Main parameter class that contains the various other sub-parameter classes.
double compute_conservation(std::shared_ptr< PHiLiP::DGBase< dim, double > > &dg, const double poly_degree) const
Function computes the conservation.
ODESolverParam ode_solver_param
Contains parameters for ODE solver.
double initial_time_step
Time step used in ODE solver.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
static std::shared_ptr< DGBase< dim, real, MeshType > > create_discontinuous_galerkin(const Parameters::AllParameters *const parameters_input, const unsigned int degree, const unsigned int max_degree_input, const unsigned int grid_degree_input, const std::shared_ptr< Triangulation > triangulation_input)
Creates a derived object DG, but returns it as DGBase.
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
ConvectiveNumericalFlux conv_num_flux_type
Store convective flux type.
BurgersEnergyStability(const Parameters::AllParameters *const parameters_input)
Constructor.
int run_test() const override
Run test.
dealii::ConditionalOStream pcout
ConditionalOStream.
DGBase is independent of the number of state 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.
static std::shared_ptr< ODESolverBase< dim, real, MeshType > > create_ODESolver(std::shared_ptr< DGBase< dim, real, MeshType > > dg_input)
Creates either implicit or explicit ODE solver based on parameter value(no POD basis given) ...
Projection operator corresponding to basis functions onto M-norm (L2).
static void set_initial_condition(std::shared_ptr< InitialConditionFunction< dim, nstate, double > > initial_condition_function_input, std::shared_ptr< PHiLiP::DGBase< dim, real > > dg_input, const Parameters::AllParameters *const parameters_input)
Applies the given initial condition function to the given dg object.
Base class of all the tests.
static std::shared_ptr< InitialConditionFunction< dim, nstate, real > > create_InitialConditionFunction(Parameters::AllParameters const *const param)
Construct InitialConditionFunction object from global parameter file.