[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
euler_cylinder.cpp
1 #include <stdlib.h>
2 #include <iostream>
3 
4 #include <deal.II/base/convergence_table.h>
5 
6 #include <deal.II/dofs/dof_tools.h>
7 
8 #include <deal.II/distributed/solution_transfer.h>
9 #include <deal.II/distributed/tria.h>
10 
11 #include <deal.II/grid/grid_generator.h>
12 #include <deal.II/grid/grid_refinement.h>
13 #include <deal.II/grid/grid_tools.h>
14 #include <deal.II/grid/grid_out.h>
15 #include <deal.II/grid/grid_in.h>
16 
17 #include <deal.II/numerics/vector_tools.h>
18 
19 #include <deal.II/fe/fe_values.h>
20 
21 #include <deal.II/fe/mapping_q.h>
22 
23 
24 #include "euler_cylinder.h"
25 
26 #include "physics/initial_conditions/initial_condition_function.h"
27 #include "physics/euler.h"
28 #include "physics/manufactured_solution.h"
29 #include "dg/dg_factory.hpp"
30 #include "ode_solver/ode_solver_factory.h"
31 
32 
33 namespace PHiLiP {
34 namespace Tests {
35 
36 dealii::Point<2> center(0.0,0.0);
37 const double inner_radius = 1, outer_radius = inner_radius*20;
38 
39 dealii::Point<2> warp_cylinder (const dealii::Point<2> &p)//, const double inner_radius, const double outer_radius)
40 {
41  const double rectangle_height = 1.0;
42  //const double original_radius = std::abs(p[0]);
43  const double angle = p[1]/rectangle_height * dealii::numbers::PI;
44 
45  //const double radius = std::abs(p[0]);
46 
47  const double power = 1.8;
48  const double radius = outer_radius*(inner_radius/outer_radius + pow(std::abs(p[0]), power));
49 
50  dealii::Point<2> q = p;
51  q[0] = -radius*cos(angle);
52  q[1] = radius*sin(angle);
53  return q;
54 }
55 
56 void half_cylinder(dealii::parallel::distributed::Triangulation<2> & tria,
57  const unsigned int n_cells_circle,
58  const unsigned int n_cells_radial)
59 {
60  //const double pi = dealii::numbers::PI;
61  //double inner_circumference = inner_radius*pi;
62  //double outer_circumference = outer_radius*pi;
63  //const double rectangle_height = inner_circumference;
64  dealii::Point<2> p1(-1,0.0), p2(-0.0,1.0);
65 
66  const bool colorize = true;
67 
68  std::vector<unsigned int> n_subdivisions(2);
69  n_subdivisions[0] = n_cells_radial;
70  n_subdivisions[1] = n_cells_circle;
71  dealii::GridGenerator::subdivided_hyper_rectangle (tria, n_subdivisions, p1, p2, colorize);
72 
73  dealii::GridTools::transform (&warp_cylinder, tria);
74 
75  tria.set_all_manifold_ids(0);
76  tria.set_manifold(0, dealii::SphericalManifold<2>(center));
77 
78  // Assign BC
79  for (auto cell = tria.begin_active(); cell != tria.end(); ++cell) {
80  //if (!cell->is_locally_owned()) continue;
81  for (unsigned int face=0; face<dealii::GeometryInfo<2>::faces_per_cell; ++face) {
82  if (cell->face(face)->at_boundary()) {
83  unsigned int current_id = cell->face(face)->boundary_id();
84  if (current_id == 0) {
85  cell->face(face)->set_boundary_id (1004); // x_left, Farfield
86  } else if (current_id == 1) {
87  cell->face(face)->set_boundary_id (1001); // x_right, Symmetry/Wall
88  } else if (current_id == 2) {
89  cell->face(face)->set_boundary_id (1001); // y_bottom, Symmetry/Wall
90  } else if (current_id == 3) {
91  cell->face(face)->set_boundary_id (1001); // y_top, Wall
92  } else {
93  std::abort();
94  }
95  }
96  }
97  }
98 }
99 
100 template <int dim, int nstate>
102  :
103  TestsBase::TestsBase(parameters_input)
104 {}
105 
106 template<int dim, int nstate>
108 ::run_test () const
109 {
110  pcout << " Running Euler cylinder entropy convergence. " << std::endl;
111  static_assert(dim==2);
113  using GridEnum = ManParam::GridEnum;
115 
116  Assert(dim == param.dimension, dealii::ExcDimensionMismatch(dim, param.dimension));
117  Assert(param.pde_type == param.PartialDifferentialEquation::euler, dealii::ExcMessage("Can't run Euler case if PDE is not Euler"));
118  //if (param.pde_type == param.PartialDifferentialEquation::euler) return 1;
119 
120  ManParam manu_grid_conv_param = param.manufactured_convergence_study_param;
121 
122  const unsigned int p_start = manu_grid_conv_param.degree_start;
123  const unsigned int p_end = manu_grid_conv_param.degree_end;
124 
125  const unsigned int n_grids_input = manu_grid_conv_param.number_of_grids;
126 
127  Physics::Euler<dim,nstate,double> euler_physics_double
129  &param,
130  param.euler_param.ref_length,
131  param.euler_param.gamma_gas,
132  param.euler_param.mach_inf,
133  param.euler_param.angle_of_attack,
134  param.euler_param.side_slip_angle);
135  FreeStreamInitialConditions<dim,nstate,double> initial_conditions(euler_physics_double);
136 
137  std::vector<int> fail_conv_poly;
138  std::vector<double> fail_conv_slop;
139  std::vector<dealii::ConvergenceTable> convergence_table_vector;
140 
141  for (unsigned int poly_degree = p_start; poly_degree <= p_end; ++poly_degree) {
142 
143  // p0 tends to require a finer grid to reach asymptotic region
144  unsigned int n_grids = n_grids_input;
145  if (poly_degree <= 1) n_grids = n_grids_input;
146 
147  std::vector<double> entropy_error(n_grids);
148  std::vector<double> grid_size(n_grids);
149 
150  const std::vector<int> n_1d_cells = get_number_1d_cells(n_grids);
151 
152  dealii::ConvergenceTable convergence_table;
153 
154  // Generate grid and mapping
155  using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
156  std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
157  MPI_COMM_WORLD,
158  typename dealii::Triangulation<dim>::MeshSmoothing(
159  dealii::Triangulation<dim>::smoothing_on_refinement |
160  dealii::Triangulation<dim>::smoothing_on_coarsening));
161 
162  const unsigned int n_cells_circle = n_1d_cells[0];
163  const unsigned int n_cells_radial = 1.5*n_cells_circle;
164  half_cylinder(*grid, n_cells_circle, n_cells_radial);
165 
166  // Create DG object
167  const int solution_degree = poly_degree;
168  //const int grid_degree = std::max(2,solution_degree);
169  const int grid_degree = solution_degree+1;
170  std::shared_ptr < DGBase<dim, double> > dg = DGFactory<dim,double>::create_discontinuous_galerkin(&param, solution_degree, solution_degree, grid_degree, grid);
171 
172  dg->allocate_system ();
173  // Initialize coarse grid solution with free-stream
174  dealii::VectorTools::interpolate(dg->dof_handler, initial_conditions, dg->solution);
175 
176  // Create ODE solver and ramp up the solution from p0
177  std::shared_ptr<ODE::ODESolverBase<dim, double>> ode_solver = ODE::ODESolverFactory<dim, double>::create_ODESolver(dg);
178  ode_solver->initialize_steady_polynomial_ramping(poly_degree);
179 
180  dealii::Vector<float> estimated_error_per_cell(grid->n_active_cells());
181  for (unsigned int igrid=0; igrid<n_grids; ++igrid) {
182 
183  // Interpolate solution from previous grid
184  if (igrid>0) {
185  dealii::LinearAlgebra::distributed::Vector<double> old_solution(dg->solution);
186  dealii::parallel::distributed::SolutionTransfer<dim, dealii::LinearAlgebra::distributed::Vector<double>, dealii::DoFHandler<dim>> solution_transfer(dg->dof_handler);
187  solution_transfer.prepare_for_coarsening_and_refinement(old_solution);
188  dg->high_order_grid->prepare_for_coarsening_and_refinement();
189 
190  grid->refine_global (1);
191  //dealii::GridRefinement::refine_and_coarsen_fixed_number(*grid,
192  // estimated_error_per_cell,
193  // 0.3,
194  // 0.03);
195  //grid->execute_coarsening_and_refinement();
196  dg->high_order_grid->execute_coarsening_and_refinement();
197  dg->allocate_system ();
198  dg->solution.zero_out_ghosts();
199  solution_transfer.interpolate(dg->solution);
200  dg->solution.update_ghost_values();
201  }
202 
203  // std::string filename = "grid_cylinder-" + dealii::Utilities::int_to_string(igrid, 1) + ".eps";
204  // std::ofstream out (filename);
205  // dealii::GridOut grid_out;
206  // grid_out.write_eps (*grid, out);
207  // pcout << " Grid #" << igrid+1 << " . Number of cells: " << grid->n_global_active_cells() << std::endl;
208  // pcout << " written to " << filename << std::endl << std::endl;
209 
210 
211  const unsigned int n_global_active_cells = grid->n_global_active_cells();
212  const unsigned int n_dofs = dg->dof_handler.n_dofs();
213  pcout << "Dimension: " << dim << "\t Polynomial degree p: " << poly_degree << std::endl
214  << "Grid number: " << igrid+1 << "/" << n_grids
215  << ". Number of active cells: " << n_global_active_cells
216  << ". Number of degrees of freedom: " << n_dofs
217  << std::endl;
218 
219  //ode_solver->initialize_steady_polynomial_ramping (poly_degree);
220  // Solve the steady state problem
221  ode_solver->steady_state();
222 
223  const auto mapping = (*(dg->high_order_grid->mapping_fe_field));
224  dealii::hp::MappingCollection<dim> mapping_collection(mapping);
225  dealii::hp::FEValues<dim,dim> fe_values_collection_volume (mapping_collection, dg->fe_collection, dg->volume_quadrature_collection, dealii::update_values | dealii::update_JxW_values);
226  // Overintegrate the error to make sure there is not integration error in the error estimate
227  //int overintegrate = 0;
228  //dealii::QGauss<dim> quad_extra(dg->max_degree+1+overintegrate);
229  //dealii::FEValues<dim,dim> fe_values_extra(*(dg->high_order_grid->mapping_fe_field), dg->fe_collection[poly_degree], quad_extra,
230  // dealii::update_values | dealii::update_JxW_values | dealii::update_quadrature_points);
231  std::array<double,nstate> soln_at_q;
232 
233  double l2error = 0;
234  double area = 0;
235  const double exact_area = (std::pow(outer_radius+inner_radius, 2.0) - std::pow(inner_radius,2.0))*dealii::numbers::PI / 2.0;
236 
237 
238  const double entropy_inf = euler_physics_double.entropy_inf;
239 
240  estimated_error_per_cell.reinit(grid->n_active_cells());
241  // Integrate solution error and output error
242  for (auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
243 
244  if (!cell->is_locally_owned()) continue;
245 
246  const int i_fele = cell->active_fe_index();
247  const int i_quad = i_fele;
248  const int i_mapp = 0;
249 
250  fe_values_collection_volume.reinit (cell, i_quad, i_mapp, i_fele);
251  const dealii::FEValues<dim,dim> &fe_values_volume = fe_values_collection_volume.get_present_fe_values();
252 
253  //fe_values_volume.reinit (cell);
254  std::vector<dealii::types::global_dof_index> dofs_indices (fe_values_volume.dofs_per_cell);
255  cell->get_dof_indices (dofs_indices);
256 
257  double cell_l2error = 0;
258  for (unsigned int iquad=0; iquad<fe_values_volume.n_quadrature_points; ++iquad) {
259 
260  std::fill(soln_at_q.begin(), soln_at_q.end(), 0);
261  for (unsigned int idof=0; idof<fe_values_volume.dofs_per_cell; ++idof) {
262  const unsigned int istate = fe_values_volume.get_fe().system_to_component_index(idof).first;
263  soln_at_q[istate] += dg->solution[dofs_indices[idof]] * fe_values_volume.shape_value_component(idof, iquad, istate);
264  }
265  const double entropy = euler_physics_double.compute_entropy_measure(soln_at_q);
266 
267  const double uexact = entropy_inf;
268  cell_l2error += pow(entropy - uexact, 2) * fe_values_volume.JxW(iquad);
269  estimated_error_per_cell[cell->active_cell_index()] = cell_l2error;
270  l2error += cell_l2error;
271 
272  area += fe_values_volume.JxW(iquad);
273  }
274  }
275  const double l2error_mpi_sum = std::sqrt(dealii::Utilities::MPI::sum(l2error, mpi_communicator));
276  const double area_mpi_sum = dealii::Utilities::MPI::sum(area, mpi_communicator);
277 
278  // Convergence table
279  const double dx = 1.0/pow(n_dofs,(1.0/dim));
280  grid_size[igrid] = dx;
281  entropy_error[igrid] = l2error_mpi_sum;
282 
283  convergence_table.add_value("p", poly_degree);
284  convergence_table.add_value("cells", n_global_active_cells);
285  convergence_table.add_value("DoFs", n_dofs);
286  convergence_table.add_value("dx", dx);
287  convergence_table.add_value("L2_entropy_error", l2error_mpi_sum);
288  convergence_table.add_value("area_error", std::abs(area_mpi_sum-exact_area));
289 
290 
291  pcout << " Grid size h: " << dx
292  << " L2-entropy_error: " << l2error_mpi_sum
293  << " Residual: " << ode_solver->residual_norm
294  << std::endl;
295 
296  if (igrid > 0) {
297  const double slope_soln_err = log(entropy_error[igrid]/entropy_error[igrid-1])
298  / log(grid_size[igrid]/grid_size[igrid-1]);
299  pcout << "From grid " << igrid-1
300  << " to grid " << igrid
301  << " dimension: " << dim
302  << " polynomial degree p: " << poly_degree
303  << std::endl
304  << " entropy_error1 " << entropy_error[igrid-1]
305  << " entropy_error2 " << entropy_error[igrid]
306  << " slope " << slope_soln_err
307  << std::endl;
308  }
309 
310  //output_results (igrid);
311  }
312  pcout << " ********************************************" << std::endl
313  << " Convergence rates for p = " << poly_degree << std::endl
314  << " ********************************************" << std::endl;
315  convergence_table.evaluate_convergence_rates("L2_entropy_error", "cells", dealii::ConvergenceTable::reduction_rate_log2, dim);
316  convergence_table.evaluate_convergence_rates("area_error", "cells", dealii::ConvergenceTable::reduction_rate_log2, dim);
317  convergence_table.set_scientific("dx", true);
318  convergence_table.set_scientific("L2_entropy_error", true);
319  convergence_table.set_scientific("area_error", true);
320  if (pcout.is_active()) convergence_table.write_text(pcout.get_stream());
321 
322  convergence_table_vector.push_back(convergence_table);
323 
324  const double expected_slope = poly_degree+1;
325 
326  const double last_slope = log(entropy_error[n_grids-1]/entropy_error[n_grids-2])
327  / log(grid_size[n_grids-1]/grid_size[n_grids-2]);
328  double before_last_slope = last_slope;
329  if ( n_grids > 2 ) {
330  before_last_slope = log(entropy_error[n_grids-2]/entropy_error[n_grids-3])
331  / log(grid_size[n_grids-2]/grid_size[n_grids-3]);
332  }
333  const double slope_avg = 0.5*(before_last_slope+last_slope);
334  const double slope_diff = slope_avg-expected_slope;
335 
336  double slope_deficit_tolerance = -std::abs(manu_grid_conv_param.slope_deficit_tolerance);
337  if(poly_degree == 0) slope_deficit_tolerance *= 2; // Otherwise, grid sizes need to be much bigger for p=0
338 
339  if (slope_diff < slope_deficit_tolerance) {
340  pcout << std::endl
341  << "Convergence order not achieved. Average last 2 slopes of "
342  << slope_avg << " instead of expected "
343  << expected_slope << " within a tolerance of "
344  << slope_deficit_tolerance
345  << std::endl;
346  // p=0 just requires too many meshes to get into the asymptotic region.
347  if(poly_degree!=0) fail_conv_poly.push_back(poly_degree);
348  if(poly_degree!=0) fail_conv_slop.push_back(slope_avg);
349  }
350 
351  }
352  pcout << std::endl << std::endl << std::endl << std::endl;
353  pcout << " ********************************************" << std::endl;
354  pcout << " Convergence summary" << std::endl;
355  pcout << " ********************************************" << std::endl;
356  for (auto conv = convergence_table_vector.begin(); conv!=convergence_table_vector.end(); conv++) {
357  if (pcout.is_active()) conv->write_text(pcout.get_stream());
358  pcout << " ********************************************" << std::endl;
359  }
360  int n_fail_poly = fail_conv_poly.size();
361  if (n_fail_poly > 0) {
362  for (int ifail=0; ifail < n_fail_poly; ++ifail) {
363  const double expected_slope = fail_conv_poly[ifail]+1;
364  const double slope_deficit_tolerance = -0.1;
365  pcout << std::endl
366  << "Convergence order not achieved for polynomial p = "
367  << fail_conv_poly[ifail]
368  << ". Slope of "
369  << fail_conv_slop[ifail] << " instead of expected "
370  << expected_slope << " within a tolerance of "
371  << slope_deficit_tolerance
372  << std::endl;
373  }
374  }
375  return n_fail_poly;
376 }
377 
378 #if PHILIP_DIM==2
380 #endif
381 
382 } // Tests namespace
383 } // PHiLiP namespace
384 
385 
unsigned int degree_start
First polynomial degree to start the loop. If diffusion, must be at least 1.
PartialDifferentialEquation pde_type
Store the PDE type to be solved.
EulerCylinder(const Parameters::AllParameters *const parameters_input)
Constructor.
const double mach_inf
Farfield Mach number.
Definition: euler.h:113
Function used to evaluate farfield conservative solution.
const double angle_of_attack
Angle of attack.
Definition: euler.h:118
Parameters related to the manufactured convergence study.
const MPI_Comm mpi_communicator
MPI communicator.
Definition: tests.h:39
Files for the baseline physics.
Definition: ADTypes.hpp:10
const double side_slip_angle
Sideslip angle.
Definition: euler.h:122
unsigned int dimension
Number of dimensions. Note that it has to match the executable PHiLiP_xD.
Main parameter class that contains the various other sub-parameter classes.
ManufacturedConvergenceStudyParam manufactured_convergence_study_param
Contains parameters for manufactured convergence study.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
Definition: tests.h:20
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.
Definition: dg_factory.cpp:10
const double entropy_inf
Entropy measure at infinity.
Definition: euler.h:127
int run_test() const override
Grid convergence on Euler Gaussian Bump.
Performs grid convergence for various polynomial degrees.
dealii::ConditionalOStream pcout
ConditionalOStream.
Definition: tests.h:45
std::vector< int > get_number_1d_cells(const int ngrids) const
Evaluates the number of cells to generate the grids for 1D grid based on input file.
Definition: tests.cpp:71
const double ref_length
Reference length.
Definition: euler.h:104
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) ...
Base class of all the tests.
Definition: tests.h:17
real compute_entropy_measure(const std::array< real, nstate > &conservative_soln) const
Evaluate entropy from conservative variables.
Definition: euler.cpp:288