1 #include "advection_explicit_periodic.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/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> 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" 35 template <
int dim,
int nstate>
40 template<
int dim,
int nstate>
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);
48 dg->global_mass_matrix.vmult( mass_matrix_times_solution, dg->solution);
51 energy = dg->solution * mass_matrix_times_solution;
56 template<
int dim,
int nstate>
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);
65 dg->global_mass_matrix.vmult( mass_matrix_times_solution, dg->solution);
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);
75 std::vector<double> ones_hat(n_dofs_cell);
78 vol_projection.
build_1D_volume_operator(dg->oneD_fe_collection[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
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];
92 conservation = ones_hat_global * mass_matrix_times_solution;
97 template <
int dim,
int nstate>
101 printf(
"starting test\n");
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 >;
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;
118 dealii::ConvergenceTable convergence_table;
119 const unsigned int igrid_start = 3;
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));
129 using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
130 std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
132 typename dealii::Triangulation<dim>::MeshSmoothing(
133 dealii::Triangulation<dim>::smoothing_on_refinement |
134 dealii::Triangulation<dim>::smoothing_on_coarsening));
138 PHiLiP::Grids::nonsymmetric_curved_grid<dim,Triangulation>(*grid, igrid);
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));
147 std::cout <<
"cells " <<n_global_active_cells2 << std::endl;
151 dg->allocate_system (
false,
false,
false);
153 std::cout <<
"Implement initial conditions" << std::endl;
155 std::shared_ptr< InitialConditionFunction<dim,nstate,double> > initial_condition_function =
161 double finalTime = 2.0;
162 if constexpr(dim==3) finalTime = 0.1;
169 ode_solver->current_iteration = 0;
172 ode_solver->advance_solution_time(dt/10.0);
180 std::ofstream myfile (
"energy_plot.gpl" , std::ios::trunc);
182 ode_solver->current_iteration = 0;
185 for (
int i = 0; i < std::ceil(finalTime/dt); ++ i){
187 ode_solver->advance_solution_time(dt);
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)
196 this->
pcout <<
" Energy was not monotonically decreasing" << std::endl;
199 if ( (current_energy*initial_energy - initial_energy >= 1.0e-12) && (all_parameters_new.
conv_num_flux_type == Parameters::AllParameters::ConvectiveNumericalFlux::central_flux))
201 this->
pcout <<
" Energy was not conserved" << std::endl;
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)
214 this->
pcout <<
"Not conserved" << std::endl;
230 ode_solver->current_iteration = 0;
233 ode_solver->advance_solution_time(finalTime);
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
241 <<
"Grid number: " << igrid+1 <<
"/" << n_grids
242 <<
". Number of active cells: " << n_global_active_cells
243 <<
". Number of degrees of freedom: " << n_dofs
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;
254 double l2error = 0.0;
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;
261 for (
auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
262 if (!cell->is_locally_owned())
continue;
264 fe_values_extra.reinit (cell);
265 cell->get_dof_indices (dofs_indices);
267 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
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);
275 for (
int istate=0; istate<nstate; ++istate) {
276 const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
278 for(
int idim=0; idim<dim; idim++){
280 uexact *= sin(2.0*pi*(qpoint[idim]-adv_speeds[idim]*finalTime));
283 uexact *= sin(pi*(qpoint[idim]- adv_speeds[idim]*finalTime));
286 l2error += pow(soln_at_q[istate] - uexact, 2) * fe_values_extra.JxW(iquad);
288 double inf_temp = std::abs(soln_at_q[istate]-uexact);
290 if(inf_temp > linf_error){
291 linf_error = inf_temp;
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));
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;
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);
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
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]);
329 this->
pcout <<
"From grid " << igrid-1
330 <<
" to grid " << igrid
331 <<
" dimension: " << dim
332 <<
" polynomial degree p: " << poly_degree
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
342 if(igrid == n_grids-1){
343 if(std::abs(slope_soln_err_inf-(poly_degree+1))>0.1)
348 this->
pcout <<
" ********************************************" 350 <<
" Convergence rates for p = " << poly_degree
352 <<
" ********************************************" 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());
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.
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.
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.
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.
int run_test() const override
Run test.
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
ConvectiveNumericalFlux conv_num_flux_type
Store convective flux type.
AdvectionPeriodic(const Parameters::AllParameters *const parameters_input)
Constructor.
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.
double compute_conservation(std::shared_ptr< PHiLiP::DGBase< dim, double > > &dg, const double poly_degree) const
Function computes the conservation.