[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
convection_diffusion_explicit_periodic.cpp
1 #include "convection_diffusion_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/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>
16 
17 #include <fstream>
18 
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"
27 
28 namespace PHiLiP {
29 namespace Tests {
30 
31 template <int dim, int nstate>
33  : TestsBase::TestsBase(parameters_input)
34 {}
35 
36 template<int dim, int nstate>
38 {
39  double energy = 0.0;
40  dg->assemble_residual();
41  energy = dg->solution * dg->right_hand_side;
42 
43  //diffusion contribution
44  const double diff_coeff = dg->all_parameters->manufactured_convergence_study_param.manufactured_solution_param.diffusion_coefficient;
45  const dealii::Tensor<2,3,double> diff_tensor = Parameters::ManufacturedSolutionParam::get_default_diffusion_tensor();
46  for(int idim=0; idim<dim; idim++){
47  for(int jdim=0; jdim<dim; jdim++){
48  double temp_energy = dg->auxiliary_solution[jdim] * dg->auxiliary_right_hand_side[idim] * diff_tensor[idim][jdim];
49  energy += diff_coeff * temp_energy;
50  }
51  }
52 
53  return energy;
54 }
55 
56 template<int dim, int nstate>
57 double ConvectionDiffusionPeriodic<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);
70  for(unsigned int iquad=0; iquad<n_quad_pts; iquad++){
71  ones[iquad] = 1.0;
72  }
73  // Projected vector of ones. That is, the interpolation of ones_hat to the volume nodes is 1.
74  std::vector<double> ones_hat(n_dofs_cell);
75  // We have to project the vector of ones because the mass matrix has an interpolation from solution nodes built into it.
76  OPERATOR::vol_projection_operator<dim,2*dim,double> vol_projection(dg->nstate, poly_degree, dg->max_grid_degree);
77  vol_projection.build_1D_volume_operator(dg->oneD_fe_collection[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
78  vol_projection.matrix_vector_mult_1D(ones, ones_hat, vol_projection.oneD_vol_operator);
79 
80  dealii::LinearAlgebra::distributed::Vector<double> ones_hat_global(dg->right_hand_side);
81  std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs_cell);
82  for (auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
83  if (!cell->is_locally_owned()) continue;
84  cell->get_dof_indices (dofs_indices);
85  for(unsigned int idof=0;idof<n_dofs_cell; idof++){
86  ones_hat_global[dofs_indices[idof]] = ones_hat[idof];
87  }
88  }
89 
90  conservation = ones_hat_global * mass_matrix_times_solution;
91 
92  // diffussive contribution to conservation
94  // const double diff_coeff = 1.0;
95  const dealii::Tensor<2,3,double> diff_tensor = Parameters::ManufacturedSolutionParam::get_default_diffusion_tensor();
96  for(int idim=0; idim<dim; idim++){
97  dealii::LinearAlgebra::distributed::Vector<double> mass_matrix_times_auxiliary_variable(dg->right_hand_side);
98  if(dg->all_parameters->use_inverse_mass_on_the_fly)
99  dg->apply_global_mass_matrix(dg->auxiliary_solution[idim], mass_matrix_times_auxiliary_variable,true);
100  else
101  dg->global_mass_matrix_auxiliary.vmult( mass_matrix_times_auxiliary_variable, dg->auxiliary_solution[idim]);
102  for(int jdim=0; jdim<dim; jdim++){
103  double temp_cons = ones_hat_global * mass_matrix_times_auxiliary_variable * diff_tensor[idim][jdim];
104  conservation += diff_coeff * temp_cons;
105  }
106  // double temp_conservation = ones_hat_global * mass_matrix_times_auxiliary_variable;
107  // conservation += diff_coeff * temp_conservation;
108  }
109 
110  return conservation;
111 }
112 
113 
114 template <int dim, int nstate>
116 {
117  this->pcout << " Running Convection Diffusion Periodicity test. " << std::endl;
118 
119  PHiLiP::Parameters::AllParameters all_parameters_new = *all_parameters;
120  double left = 0.0;
121  double right = 2.0;
122  const unsigned int n_grids = (all_parameters_new.use_energy) ? 4 : 5;
123  std::vector<double> grid_size(n_grids);
124  std::vector<double> soln_error(n_grids);
125  unsigned int poly_degree = 3;
126  dealii::ConvergenceTable convergence_table;
127  const unsigned int igrid_start = (all_parameters_new.use_energy) ? 3 : 3;
128  const unsigned int grid_degree = 1;
129 
130  for(unsigned int igrid = igrid_start; igrid<n_grids; igrid++) {
131 
132 #if PHILIP_DIM==1 // dealii::parallel::distributed::Triangulation<dim> does not work for 1D
133  using Triangulation = dealii::Triangulation<dim>;
134  std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
135  typename dealii::Triangulation<dim>::MeshSmoothing(
136  dealii::Triangulation<dim>::smoothing_on_refinement |
137  dealii::Triangulation<dim>::smoothing_on_coarsening));
138 #else
139  using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
140  std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
141  MPI_COMM_WORLD,
142  typename dealii::Triangulation<dim>::MeshSmoothing(
143  dealii::Triangulation<dim>::smoothing_on_refinement |
144  dealii::Triangulation<dim>::smoothing_on_coarsening));
145 #endif
146  // straight grid
147  dealii::GridGenerator::hyper_cube(*grid, left, right, true);
148 #if PHILIP_DIM==1
149  std::vector<dealii::GridTools::PeriodicFacePair<typename dealii::Triangulation<PHILIP_DIM>::cell_iterator> > matched_pairs;
150  dealii::GridTools::collect_periodic_faces(*grid,0,1,0,matched_pairs);
151  grid->add_periodicity(matched_pairs);
152 #else
153  std::vector<dealii::GridTools::PeriodicFacePair<typename dealii::parallel::distributed::Triangulation<PHILIP_DIM>::cell_iterator> > matched_pairs;
154  dealii::GridTools::collect_periodic_faces(*grid,0,1,0,matched_pairs);
155  if(dim>=2) dealii::GridTools::collect_periodic_faces(*grid,2,3,1,matched_pairs);
156  if(dim>=3) dealii::GridTools::collect_periodic_faces(*grid,4,5,2,matched_pairs);
157  grid->add_periodicity(matched_pairs);
158 #endif
159  grid->refine_global(igrid);
160  this->pcout << "Grid generated and refined" << std::endl;
161  // CFL number
162  const unsigned int n_global_active_cells2 = grid->n_global_active_cells();
163  double n_dofs_cfl = pow(n_global_active_cells2,dim) * pow(poly_degree+1.0, dim);
164  double delta_x = (right-left)/pow(n_dofs_cfl,(1.0/dim));
166  //const double diff_coeff2 = dg->all_parameters->manufactured_convergence_study_param.manufactured_solution_param.diffusion_coefficient;
167  const dealii::Tensor<2,3,double> diff_tensor2 = Parameters::ManufacturedSolutionParam::get_default_diffusion_tensor();
168  double max_diff_tens=0.0;
169  for(int i=0; i<dim; i++){
170  for(int j=0; j<dim; j++){
171  if(std::abs(diff_tensor2[i][j]) > max_diff_tens) max_diff_tens = std::abs(diff_tensor2[i][j]);
172  }
173  }
174  all_parameters_new.ode_solver_param.initial_time_step = 0.5*pow(delta_x,2)/diff_coeff2;
175  all_parameters_new.ode_solver_param.initial_time_step = 0.05*pow(delta_x,2)/diff_coeff2 / max_diff_tens;
176 
177  // allocate dg
178  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);
179  this->pcout << "dg created" <<std::endl;
180  dg->allocate_system (false,false,false);
181 
182  this->pcout << "Setting up Initial Condition" << std::endl;
183  // Create initial condition function
184  std::shared_ptr< InitialConditionFunction<dim,nstate,double> > initial_condition_function =
186  SetInitialCondition<dim,nstate,double>::set_initial_condition(initial_condition_function, dg, &all_parameters_new);
187 
188  // Create ODE solver using the factory and providing the DG object
189  std::shared_ptr<ODE::ODESolverBase<dim, double>> ode_solver = ODE::ODESolverFactory<dim, double>::create_ODESolver(dg);
190 
191  double finalTime = 2.0;
192 
193  if(all_parameters_new.use_energy == true){//for split form get energy
194  finalTime = 0.01;
195 
196  double dt = all_parameters_new.ode_solver_param.initial_time_step;
197 
198  // need to call ode_solver before calculating energy because mass matrix isn't allocated yet.
199 
200  ode_solver->current_iteration = 0;
201  ode_solver->advance_solution_time(0.000001);
202 
203  double initial_energy = compute_energy_derivative_norm(dg);
204  double initial_conservation = compute_conservation(dg, poly_degree);
205 
206  // currently the only way to calculate energy at each time-step is to advance solution by dt instead of finaltime
207  // this causes some issues with outputs (only one file is output, which is overwritten at each time step)
208  // also the ode solver output doesn't make sense (says "iteration 1 out of 1")
209  // but it works. I'll keep it for now and need to modify the output functions later to account for this.
210  std::ofstream myfile ("energy_plot_cons_diff.gpl" , std::ios::trunc);
211 
212  ode_solver->current_iteration = 0;
213  for(int i = 0; i < std::ceil(finalTime/dt); ++ i) {
214  ode_solver->advance_solution_time(dt);
215  double current_energy = compute_energy_derivative_norm(dg);
216  std::cout << std::setprecision(16) << std::fixed;
217  this->pcout << "Energy at time " << i * dt << " is " << current_energy << std::endl;
218  myfile << i * dt << " " << std::fixed << std::setprecision(16) << current_energy << std::endl;
219 
220  if (current_energy>1e-12)/*if (current_energy*initial_energy - initial_energy >= 10000.0)*/
221  {
222  this->pcout<<"Energy Fail Not montonically decrease with current "<< current_energy<<" vs initial "<<initial_energy<<std::endl;
223  return 1;
224  break;
225  }
226 
227  if (std::abs(current_energy - initial_energy >= 1.0e-12) &&
228  (all_parameters_new.diss_num_flux_type == Parameters::AllParameters::DissipativeNumericalFlux::central_visc_flux))
229  {
230  this->pcout<<"Energy Not conserved with current "<< current_energy<<" vs initial "<<initial_energy<<std::endl;
231  return 1;
232  break;
233  }
234 
235  double current_conservation = compute_conservation(dg, poly_degree);
236  current_conservation /=initial_conservation;
237 
238  std::cout << std::setprecision(16) << std::fixed;
239  this->pcout << "Normalized Conservation at time " << i * dt << " is " << current_conservation<< std::endl;
240  myfile << i * dt << " " << std::fixed << std::setprecision(16) << current_conservation << std::endl;
241  if (current_conservation*initial_conservation - initial_conservation >= 1e-12)
242  //if (current_energy - initial_energy >= 10.00)
243  {
244  this->pcout << "Not conserved" << std::endl;
245  return 1;
246  break;
247  }
248  }
249  myfile.close();
250  }//end of energy
251  else {//do OOA
252 
253  finalTime = 1e-3;
254 
255  ode_solver->current_iteration = 0;
256 
257  ode_solver->advance_solution_time(finalTime);
258  const unsigned int n_global_active_cells = grid->n_global_active_cells();
259  const unsigned int n_dofs = dg->dof_handler.n_dofs();
260  this->pcout << "Dimension: " << dim
261  << "\t Polynomial degree p: " << poly_degree
262  << std::endl
263  << "Grid number: " << igrid+1 << "/" << n_grids
264  << ". Number of active cells: " << n_global_active_cells
265  << ". Number of degrees of freedom: " << n_dofs
266  << std::endl;
267 
268  // Overintegrate the error to make sure there is not integration error in the error estimate
269  int overintegrate = 10;
270  dealii::QGauss<dim> quad_extra(poly_degree+1+overintegrate);
271  dealii::FEValues<dim,dim> fe_values_extra(*(dg->high_order_grid->mapping_fe_field), dg->fe_collection[poly_degree], quad_extra,
272  dealii::update_values | dealii::update_JxW_values | dealii::update_quadrature_points);
273  const unsigned int n_quad_pts = fe_values_extra.n_quadrature_points;
274  std::array<double,nstate> soln_at_q;
275 
276  double l2error = 0.0;
277 
278  // Integrate solution error and output error
279  const double pi = atan(1)*4.0;
281  std::vector<dealii::types::global_dof_index> dofs_indices (fe_values_extra.dofs_per_cell);
282  for (auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
283 
284  if (!cell->is_locally_owned()) continue;
285 
286  fe_values_extra.reinit (cell);
287  cell->get_dof_indices (dofs_indices);
288 
289  for (unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
290 
291  std::fill(soln_at_q.begin(), soln_at_q.end(), 0.0);
292  for (unsigned int idof=0; idof<fe_values_extra.dofs_per_cell; ++idof) {
293  const unsigned int istate = fe_values_extra.get_fe().system_to_component_index(idof).first;
294  soln_at_q[istate] += dg->solution[dofs_indices[idof]] * fe_values_extra.shape_value_component(idof, iquad, istate);
295  }
296 
297  for (int istate=0; istate<nstate; ++istate) {
298  const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
299  double uexact=1.0;
300  for(int idim=0; idim<dim; idim++){
301  uexact *= sin(pi*(qpoint[idim]));//for grid 1-3
302  }
303  uexact *= exp(- diff_coeff * finalTime);
304  l2error += pow(soln_at_q[istate] - uexact, 2) * fe_values_extra.JxW(iquad);
305  }
306  }
307  }
308 #if PHILIP_DIM==1
309  const double l2error_mpi_sum = sqrt(l2error);
310 #else
311  const double l2error_mpi_sum = std::sqrt(dealii::Utilities::MPI::sum(l2error, this->mpi_communicator));
312 #endif
313 
314  // Convergence table
315  const double dx = 1.0/pow(n_dofs,(1.0/dim));
316  grid_size[igrid] = dx;
317  soln_error[igrid] = l2error_mpi_sum;
318  convergence_table.add_value("p", poly_degree);
319  convergence_table.add_value("cells", n_global_active_cells);
320  convergence_table.add_value("DoFs", n_dofs);
321  convergence_table.add_value("dx", dx);
322  convergence_table.add_value("soln_L2_error", l2error_mpi_sum);
323 
324  this->pcout << " Grid size h: " << dx
325  << " L2-soln_error: " << l2error_mpi_sum
326  << " Residual: " << ode_solver->residual_norm
327  << std::endl;
328 
329  if (igrid > igrid_start) {
330  const double slope_soln_err = log(soln_error[igrid]/soln_error[igrid-1])
331  / log(grid_size[igrid]/grid_size[igrid-1]);
332  // const double slope_output_err = log(output_error[igrid]/output_error[igrid-1])
333  // / log(grid_size[igrid]/grid_size[igrid-1]);
334  this->pcout << "From grid " << igrid
335  << " to grid " << igrid+1
336  << " dimension: " << dim
337  << " polynomial degree p: " << poly_degree
338  << std::endl
339  << " solution_error1 " << soln_error[igrid-1]
340  << " solution_error2 " << soln_error[igrid]
341  << " slope " << slope_soln_err
342  << std::endl;
343  if(igrid == n_grids-1){
344  if(std::abs(slope_soln_err-(poly_degree+1))>0.1){
345  return 1;
346  }
347  }
348  }
349 
350  this->pcout << " ********************************************"
351  << std::endl
352  << " Convergence rates for p = " << poly_degree
353  << std::endl
354  << " ********************************************"
355  << std::endl;
356  convergence_table.evaluate_convergence_rates("soln_L2_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  if (this->pcout.is_active()) convergence_table.write_text(this->pcout.get_stream());
360  }//end of OOA
361  }//end of grid loop
362 
363  return 0;
364 }
365 
367 
368 } // Tests namespace
369 } // PHiLiP namespace
double compute_conservation(std::shared_ptr< PHiLiP::DGBase< dim, double > > &dg, const double poly_degree) const
Function computes the conservation.
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
double compute_energy_derivative_norm(std::shared_ptr< PHiLiP::DGBase< dim, double > > &dg) const
Function computes the energy bound.
bool use_energy
Flag to use an energy monotonicity test.
Files for the baseline physics.
Definition: ADTypes.hpp:10
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 dealii::Tensor< 2, 3, double > get_default_diffusion_tensor()
gets the default diffusion tensor
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
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
Definition: operators.h:355
Convection Diffusion periodic unsteady test (currently only diffusion)
static double get_default_diffusion_coefficient()
gets the default diffusion coefficient;
DissipativeNumericalFlux diss_num_flux_type
Store diffusive flux type.
dealii::ConditionalOStream pcout
ConditionalOStream.
Definition: tests.h:45
ConvectionDiffusionPeriodic(const Parameters::AllParameters *const parameters_input)
Constructor.
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.