1 #include "convection_diffusion_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/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 dg->assemble_residual();
41 energy = dg->solution * dg->right_hand_side;
44 const double diff_coeff = dg->all_parameters->manufactured_convergence_study_param.manufactured_solution_param.diffusion_coefficient;
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;
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);
70 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
74 std::vector<double> ones_hat(n_dofs_cell);
77 vol_projection.
build_1D_volume_operator(dg->oneD_fe_collection[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
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];
90 conservation = ones_hat_global * mass_matrix_times_solution;
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);
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;
114 template <
int dim,
int nstate>
117 this->
pcout <<
" Running Convection Diffusion Periodicity test. " << std::endl;
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;
130 for(
unsigned int igrid = igrid_start; igrid<n_grids; igrid++) {
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));
139 using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
140 std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
142 typename dealii::Triangulation<dim>::MeshSmoothing(
143 dealii::Triangulation<dim>::smoothing_on_refinement |
144 dealii::Triangulation<dim>::smoothing_on_coarsening));
147 dealii::GridGenerator::hyper_cube(*grid, left, right,
true);
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);
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);
159 grid->refine_global(igrid);
160 this->
pcout <<
"Grid generated and refined" << std::endl;
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));
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]);
179 this->
pcout <<
"dg created" <<std::endl;
180 dg->allocate_system (
false,
false,
false);
182 this->
pcout <<
"Setting up Initial Condition" << std::endl;
184 std::shared_ptr< InitialConditionFunction<dim,nstate,double> > initial_condition_function =
191 double finalTime = 2.0;
200 ode_solver->current_iteration = 0;
201 ode_solver->advance_solution_time(0.000001);
210 std::ofstream myfile (
"energy_plot_cons_diff.gpl" , std::ios::trunc);
212 ode_solver->current_iteration = 0;
213 for(
int i = 0; i < std::ceil(finalTime/dt); ++ i) {
214 ode_solver->advance_solution_time(dt);
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;
220 if (current_energy>1e-12)
222 this->
pcout<<
"Energy Fail Not montonically decrease with current "<< current_energy<<
" vs initial "<<initial_energy<<std::endl;
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))
230 this->
pcout<<
"Energy Not conserved with current "<< current_energy<<
" vs initial "<<initial_energy<<std::endl;
236 current_conservation /=initial_conservation;
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)
244 this->
pcout <<
"Not conserved" << std::endl;
255 ode_solver->current_iteration = 0;
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
263 <<
"Grid number: " << igrid+1 <<
"/" << n_grids
264 <<
". Number of active cells: " << n_global_active_cells
265 <<
". Number of degrees of freedom: " << n_dofs
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;
276 double l2error = 0.0;
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) {
284 if (!cell->is_locally_owned())
continue;
286 fe_values_extra.reinit (cell);
287 cell->get_dof_indices (dofs_indices);
289 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
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);
297 for (
int istate=0; istate<nstate; ++istate) {
298 const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
300 for(
int idim=0; idim<dim; idim++){
301 uexact *= sin(pi*(qpoint[idim]));
303 uexact *= exp(- diff_coeff * finalTime);
304 l2error += pow(soln_at_q[istate] - uexact, 2) * fe_values_extra.JxW(iquad);
309 const double l2error_mpi_sum = sqrt(l2error);
311 const double l2error_mpi_sum = std::sqrt(dealii::Utilities::MPI::sum(l2error, this->
mpi_communicator));
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);
324 this->
pcout <<
" Grid size h: " << dx
325 <<
" L2-soln_error: " << l2error_mpi_sum
326 <<
" Residual: " << ode_solver->residual_norm
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]);
334 this->
pcout <<
"From grid " << igrid
335 <<
" to grid " << igrid+1
336 <<
" dimension: " << dim
337 <<
" polynomial degree p: " << poly_degree
339 <<
" solution_error1 " << soln_error[igrid-1]
340 <<
" solution_error2 " << soln_error[igrid]
341 <<
" slope " << slope_soln_err
343 if(igrid == n_grids-1){
344 if(std::abs(slope_soln_err-(poly_degree+1))>0.1){
350 this->
pcout <<
" ********************************************" 352 <<
" Convergence rates for p = " << poly_degree
354 <<
" ********************************************" 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());
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.
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
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.
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 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.
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
int run_test() const override
Run test.
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.
ConvectionDiffusionPeriodic(const Parameters::AllParameters *const parameters_input)
Constructor.
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.