[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
advection_explicit_periodic.cpp
1 #include "advection_explicit_periodic.h"
2 
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/fe/mapping_q.h>
9 #include <deal.II/grid/grid_generator.h>
10 #include <deal.II/grid/grid_in.h>
11 #include <deal.II/grid/grid_out.h>
12 #include <deal.II/grid/grid_refinement.h>
13 #include <deal.II/grid/grid_tools.h>
14 #include <deal.II/grid/manifold_lib.h>
15 #include <deal.II/numerics/data_out.h>
16 #include <deal.II/numerics/solution_transfer.h>
17 #include <deal.II/numerics/vector_tools.h>
18 #include <fenv.h>
19 
20 #include <fstream>
21 
22 #include "dg/dg_base.hpp"
23 #include "dg/dg_factory.hpp"
24 #include "mesh/grids/nonsymmetric_curved_periodic_grid.hpp"
25 #include "ode_solver/ode_solver_factory.h"
26 #include "parameters/all_parameters.h"
27 #include "parameters/parameters.h"
28 #include "physics/initial_conditions/initial_condition_function.h"
29 #include "physics/initial_conditions/set_initial_condition.h"
30 #include "physics/physics.h"
31 #include "physics/physics_factory.h"
32 
33 namespace PHiLiP {
34 namespace Tests {
35 template <int dim, int nstate>
37  : TestsBase::TestsBase(parameters_input)
38 {}
39 
40 template<int dim, int nstate>
42 {
43  double energy = 0.0;
44  dealii::LinearAlgebra::distributed::Vector<double> mass_matrix_times_solution(dg->right_hand_side);
45  if(dg->all_parameters->use_inverse_mass_on_the_fly)
46  dg->apply_global_mass_matrix(dg->solution,mass_matrix_times_solution);
47  else
48  dg->global_mass_matrix.vmult( mass_matrix_times_solution, dg->solution);
49  //Since we normalize the energy later, don't bother scaling by 0.5
50  //Energy \f$ = 0.5 * \int u^2 d\Omega_m \f$
51  energy = dg->solution * mass_matrix_times_solution;
52 
53  return energy;
54 }
55 
56 template<int dim, int nstate>
57 double AdvectionPeriodic<dim, nstate>::compute_conservation(std::shared_ptr < PHiLiP::DGBase<dim, double> > &dg, const double poly_degree) const
58 {
59  //Conservation \f$ = \int 1 * u d\Omega_m \f$
60  double conservation = 0.0;
61  dealii::LinearAlgebra::distributed::Vector<double> mass_matrix_times_solution(dg->right_hand_side);
62  if(dg->all_parameters->use_inverse_mass_on_the_fly)
63  dg->apply_global_mass_matrix(dg->solution,mass_matrix_times_solution);
64  else
65  dg->global_mass_matrix.vmult( mass_matrix_times_solution, dg->solution);
66 
67  const unsigned int n_dofs_cell = dg->fe_collection[poly_degree].dofs_per_cell;
68  const unsigned int n_quad_pts = dg->volume_quadrature_collection[poly_degree].size();
69  std::vector<double> ones(n_quad_pts, 1.0);
70  // std::vector<double> ones(n_quad_pts);
71  // for(unsigned int iquad=0; iquad<n_quad_pts; iquad++){
72  // ones[iquad] = 1.0;
73  // }
74  // Projected vector of ones. That is, the interpolation of ones_hat to the volume nodes is 1.
75  std::vector<double> ones_hat(n_dofs_cell);
76  // We have to project the vector of ones because the mass matrix has an interpolation from solution nodes built into it.
77  OPERATOR::vol_projection_operator<dim,2*dim,double> vol_projection(dg->nstate, poly_degree, dg->max_grid_degree);
78  vol_projection.build_1D_volume_operator(dg->oneD_fe_collection[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
79  vol_projection.matrix_vector_mult_1D(ones, ones_hat,
80  vol_projection.oneD_vol_operator);
81 
82  dealii::LinearAlgebra::distributed::Vector<double> ones_hat_global(dg->right_hand_side);
83  std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs_cell);
84  for (auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
85  if (!cell->is_locally_owned()) continue;
86  cell->get_dof_indices (dofs_indices);
87  for(unsigned int idof=0;idof<n_dofs_cell; idof++){
88  ones_hat_global[dofs_indices[idof]] = ones_hat[idof];
89  }
90  }
91 
92  conservation = ones_hat_global * mass_matrix_times_solution;
93 
94  return conservation;
95 }
96 
97 template <int dim, int nstate>
99 {
100 
101  printf("starting test\n");
102  PHiLiP::Parameters::AllParameters all_parameters_new = *all_parameters;
103 
104  const unsigned int n_grids = (all_parameters_new.use_energy) ? 4 : 5;
105  std::vector<double> grid_size(n_grids);
106  std::vector<double> soln_error(n_grids);
107  std::vector<double> soln_error_inf(n_grids);
108  using ADtype = double;
109  using ADArray = std::array<ADtype,nstate>;
110  using ADArrayTensor1 = std::array< dealii::Tensor<1,dim,ADtype>, nstate >;
111 
112  const double left = -1.0;
113  const double right = 1.0;
114  unsigned int n_refinements = n_grids;
115  unsigned int poly_degree = 3;
116  unsigned int grid_degree = poly_degree;
117 
118  dealii::ConvergenceTable convergence_table;
119  const unsigned int igrid_start = 3;
120 
121  for(unsigned int igrid=igrid_start; igrid<n_refinements; ++igrid){
122 #if PHILIP_DIM==1 // dealii::parallel::distributed::Triangulation<dim> does not work for 1D
123  using Triangulation = dealii::Triangulation<dim>;
124  std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
125  typename dealii::Triangulation<dim>::MeshSmoothing(
126  dealii::Triangulation<dim>::smoothing_on_refinement |
127  dealii::Triangulation<dim>::smoothing_on_coarsening));
128 #else
129  using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
130  std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
131  MPI_COMM_WORLD,
132  typename dealii::Triangulation<dim>::MeshSmoothing(
133  dealii::Triangulation<dim>::smoothing_on_refinement |
134  dealii::Triangulation<dim>::smoothing_on_coarsening));
135 #endif
136 
137  //set the warped grid
138  PHiLiP::Grids::nonsymmetric_curved_grid<dim,Triangulation>(*grid, igrid);
139 
140  //CFL number
141  const unsigned int n_global_active_cells2 = grid->n_global_active_cells();
142  double n_dofs_cfl = pow(n_global_active_cells2,dim) * pow(poly_degree+1.0, dim);
143  double delta_x = (right-left)/pow(n_dofs_cfl,(1.0/dim));
144  all_parameters_new.ode_solver_param.initial_time_step = delta_x /(1.0*(2.0*poly_degree+1)) ;
145  all_parameters_new.ode_solver_param.initial_time_step = (all_parameters_new.use_energy) ? 0.05*delta_x : 0.5*delta_x;
146  std::cout << "dt " <<all_parameters_new.ode_solver_param.initial_time_step << std::endl;
147  std::cout << "cells " <<n_global_active_cells2 << std::endl;
148 
149  //Set the DG spatial sys
150  std::shared_ptr < PHiLiP::DGBase<dim, double> > dg = PHiLiP::DGFactory<dim,double>::create_discontinuous_galerkin(&all_parameters_new, poly_degree, poly_degree, grid_degree, grid);
151  dg->allocate_system (false,false,false);
152 
153  std::cout << "Implement initial conditions" << std::endl;
154  // Create initial condition function
155  std::shared_ptr< InitialConditionFunction<dim,nstate,double> > initial_condition_function =
157  SetInitialCondition<dim,nstate,double>::set_initial_condition(initial_condition_function, dg, &all_parameters_new);
158 
159  // Create ODE solver using the factory and providing the DG object
160  std::shared_ptr<PHiLiP::ODE::ODESolverBase<dim, double>> ode_solver = PHiLiP::ODE::ODESolverFactory<dim, double>::create_ODESolver(dg);
161  double finalTime = 2.0;
162  if constexpr(dim==3) finalTime = 0.1;//to speed things up locally
163 
164  // need to call ode_solver before calculating energy because mass matrix isn't allocated yet.
165 
166  if (all_parameters_new.use_energy == true){//for split form get energy
167  double dt = all_parameters_new.ode_solver_param.initial_time_step;
168 
169  ode_solver->current_iteration = 0;
170 
171  // advance by small amount, basically 0
172  ode_solver->advance_solution_time(dt/10.0);
173  // compute the initial energy and conservation
174  double initial_energy = compute_energy(dg);
175  double initial_conservation = compute_conservation(dg, poly_degree);
176 
177  // create file to write energy and conservation results
178  // outputs results as Time, Energy Newline Time, Conservation
179  // And energy and conservation values are normalized by the initial value.
180  std::ofstream myfile ("energy_plot.gpl" , std::ios::trunc);
181 
182  ode_solver->current_iteration = 0;
183 
184  // loop over time steps because needs to evaluate energy and conservation at each time step.
185  for (int i = 0; i < std::ceil(finalTime/dt); ++ i){
186 
187  ode_solver->advance_solution_time(dt);
188  //energy
189  double current_energy = compute_energy(dg);
190  current_energy /=initial_energy;
191  std::cout << std::setprecision(16) << std::fixed;
192  this->pcout << "Energy at time " << i * dt << " is " << current_energy<< std::endl;
193  myfile << i * dt << " " << std::fixed << std::setprecision(16) << current_energy << std::endl;
194  if (current_energy*initial_energy - initial_energy >= 1.00)//since normalized by initial
195  {
196  this->pcout << " Energy was not monotonically decreasing" << std::endl;
197  return 1;
198  }
199  if ( (current_energy*initial_energy - initial_energy >= 1.0e-12) && (all_parameters_new.conv_num_flux_type == Parameters::AllParameters::ConvectiveNumericalFlux::central_flux))
200  {
201  this->pcout << " Energy was not conserved" << std::endl;
202  return 1;
203  }
204 
205  // Conservation
206  double current_conservation = compute_conservation(dg, poly_degree);
207  current_conservation /=initial_conservation;
208  std::cout << std::setprecision(16) << std::fixed;
209  this->pcout << "Normalized Conservation at time " << i * dt << " is " << current_conservation<< std::endl;
210  myfile << i * dt << " " << std::fixed << std::setprecision(16) << current_conservation << std::endl;
211  if (current_conservation*initial_conservation - initial_conservation >= 10.00)
212  // if (current_energy - initial_energy >= 10.00)
213  {
214  this->pcout << "Not conserved" << std::endl;
215  return 1;
216  }
217  }
218 
219  // Close the file
220  myfile.close();
221  }
222  else{// do OOA
223  if(left==-1){
224  finalTime = 0.5;
225  }
226  if(left==0){
227  finalTime=1.0;
228  }
229 
230  ode_solver->current_iteration = 0;
231 
232  // advance solution until the final time
233  ode_solver->advance_solution_time(finalTime);
234 
235  // output results
236  const unsigned int n_global_active_cells = grid->n_global_active_cells();
237  const unsigned int n_dofs = dg->dof_handler.n_dofs();
238  this->pcout << "Dimension: " << dim
239  << "\t Polynomial degree p: " << poly_degree
240  << std::endl
241  << "Grid number: " << igrid+1 << "/" << n_grids
242  << ". Number of active cells: " << n_global_active_cells
243  << ". Number of degrees of freedom: " << n_dofs
244  << std::endl;
245 
246  // Overintegrate the error to make sure there is not integration error in the error estimate
247  int overintegrate = 10;
248  dealii::QGauss<dim> quad_extra(poly_degree+1+overintegrate);
249  dealii::FEValues<dim,dim> fe_values_extra(*(dg->high_order_grid->mapping_fe_field), dg->fe_collection[poly_degree], quad_extra,
250  dealii::update_values | dealii::update_JxW_values | dealii::update_quadrature_points);
251  const unsigned int n_quad_pts = fe_values_extra.n_quadrature_points;
252  std::array<double,nstate> soln_at_q;
253 
254  double l2error = 0.0;
255 
256  // Integrate every cell and compute L2 and Linf errors.
257  const double pi = atan(1)*4.0;
258  std::vector<dealii::types::global_dof_index> dofs_indices (fe_values_extra.dofs_per_cell);
259  double linf_error = 0.0;
260  const dealii::Tensor<1,3,double> adv_speeds = Parameters::ManufacturedSolutionParam::get_default_advection_vector();
261  for (auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
262  if (!cell->is_locally_owned()) continue;
263 
264  fe_values_extra.reinit (cell);
265  cell->get_dof_indices (dofs_indices);
266 
267  for (unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
268 
269  std::fill(soln_at_q.begin(), soln_at_q.end(), 0.0);
270  for (unsigned int idof=0; idof<fe_values_extra.dofs_per_cell; ++idof) {
271  const unsigned int istate = fe_values_extra.get_fe().system_to_component_index(idof).first;
272  soln_at_q[istate] += dg->solution[dofs_indices[idof]] * fe_values_extra.shape_value_component(idof, iquad, istate);
273  }
274 
275  for (int istate=0; istate<nstate; ++istate) {
276  const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
277  double uexact=1.0;
278  for(int idim=0; idim<dim; idim++){
279  if(left==0){
280  uexact *= sin(2.0*pi*(qpoint[idim]-adv_speeds[idim]*finalTime));//for grid 1-3
281  }
282  if(left==-1){
283  uexact *= sin(pi*(qpoint[idim]- adv_speeds[idim]*finalTime));//for grid 1-3
284  }
285  }
286  l2error += pow(soln_at_q[istate] - uexact, 2) * fe_values_extra.JxW(iquad);
287 
288  double inf_temp = std::abs(soln_at_q[istate]-uexact);
289  //store pointwise inf error
290  if(inf_temp > linf_error){
291  linf_error = inf_temp;
292  }
293  }
294  }
295 
296  }
297  const double l2error_mpi_sum = std::sqrt(dealii::Utilities::MPI::sum(l2error, this->mpi_communicator));
298  const double linferror_mpi= (dealii::Utilities::MPI::max(linf_error, this->mpi_communicator));
299 
300  // Convergence table
301  const double dx = 1.0/pow(n_dofs,(1.0/dim));
302  grid_size[igrid] = dx;
303  soln_error[igrid] = l2error_mpi_sum;
304  soln_error_inf[igrid] = linferror_mpi;
305  // output_error[igrid] = std::abs(solution_integral - exact_solution_integral);
306 
307  convergence_table.add_value("p", poly_degree);
308  convergence_table.add_value("cells", n_global_active_cells);
309  convergence_table.add_value("DoFs", n_dofs);
310  convergence_table.add_value("dx", dx);
311  convergence_table.add_value("soln_L2_error", l2error_mpi_sum);
312  convergence_table.add_value("soln_Linf_error", linferror_mpi);
313  // convergence_table.add_value("output_error", output_error[igrid]);
314 
315 
316  this->pcout << " Grid size h: " << dx
317  << " L2-soln_error: " << l2error_mpi_sum
318  << " Linf-soln_error: " << linferror_mpi
319  << " Residual: " << ode_solver->residual_norm
320  << std::endl;
321 
322  if (igrid > igrid_start) {
323  const double slope_soln_err = log(soln_error[igrid]/soln_error[igrid-1])
324  / log(grid_size[igrid]/grid_size[igrid-1]);
325  const double slope_soln_err_inf = log(soln_error_inf[igrid]/soln_error_inf[igrid-1])
326  / log(grid_size[igrid]/grid_size[igrid-1]);
327  // const double slope_output_err = log(output_error[igrid]/output_error[igrid-1])
328  // / log(grid_size[igrid]/grid_size[igrid-1]);
329  this->pcout << "From grid " << igrid-1
330  << " to grid " << igrid
331  << " dimension: " << dim
332  << " polynomial degree p: " << poly_degree
333  << std::endl
334  << " solution_error1 " << soln_error[igrid-1]
335  << " solution_error2 " << soln_error[igrid]
336  << " slope " << slope_soln_err
337  << " solution_error1_inf " << soln_error_inf[igrid-1]
338  << " solution_error2_inf " << soln_error_inf[igrid]
339  << " slope " << slope_soln_err_inf
340  << std::endl;
341 
342  if(igrid == n_grids-1){
343  if(std::abs(slope_soln_err_inf-(poly_degree+1))>0.1)
344  return 1;
345  }
346  }
347  }//end of OOA else statement
348  this->pcout << " ********************************************"
349  << std::endl
350  << " Convergence rates for p = " << poly_degree
351  << std::endl
352  << " ********************************************"
353  << std::endl;
354  convergence_table.evaluate_convergence_rates("soln_L2_error", "cells", dealii::ConvergenceTable::reduction_rate_log2, dim);
355  convergence_table.evaluate_convergence_rates("soln_Linf_error", "cells", dealii::ConvergenceTable::reduction_rate_log2, dim);
356  convergence_table.evaluate_convergence_rates("output_error", "cells", dealii::ConvergenceTable::reduction_rate_log2, dim);
357  convergence_table.set_scientific("dx", true);
358  convergence_table.set_scientific("soln_L2_error", true);
359  convergence_table.set_scientific("soln_Linf_error", true);
360  convergence_table.set_scientific("output_error", true);
361  if (this->pcout.is_active()) convergence_table.write_text(this->pcout.get_stream());
362 
363  }//end of grid loop
364 
365  return 0;//if reaches here mean passed test
366 }
367 
368 template class AdvectionPeriodic <PHILIP_DIM,1>;
369 
370 } //Tests namespace
371 } //PHiLiP namespace
double compute_energy(std::shared_ptr< PHiLiP::DGBase< dim, double > > &dg) const
Function computes the energy.
static dealii::Tensor< 1, 3, double > get_default_advection_vector()
gets the default advection vector
const MPI_Comm mpi_communicator
MPI communicator.
Definition: tests.h:39
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
Definition: operators.cpp:1732
bool use_energy
Flag to use an energy monotonicity test.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Advection periodic unsteady test.
Main parameter class that contains the various other sub-parameter classes.
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.
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
int run_test() const override
Run test.
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
Definition: operators.h:355
ConvectiveNumericalFlux conv_num_flux_type
Store convective flux type.
AdvectionPeriodic(const Parameters::AllParameters *const parameters_input)
Constructor.
dealii::ConditionalOStream pcout
ConditionalOStream.
Definition: tests.h:45
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
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.
Definition: operators.cpp:308
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).
Definition: operators.h:694
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.
Definition: tests.h:17
static std::shared_ptr< InitialConditionFunction< dim, nstate, real > > create_InitialConditionFunction(Parameters::AllParameters const *const param)
Construct InitialConditionFunction object from global parameter file.
double compute_conservation(std::shared_ptr< PHiLiP::DGBase< dim, double > > &dg, const double poly_degree) const
Function computes the conservation.