2 #include "dg/dg_factory.hpp" 3 #include "euler_split_inviscid_taylor_green_vortex.h" 4 #include "physics/initial_conditions/set_initial_condition.h" 5 #include "physics/initial_conditions/initial_condition_function.h" 6 #include "mesh/grids/nonsymmetric_curved_periodic_grid.hpp" 7 #include "mesh/grids/straight_periodic_cube.hpp" 12 template <
int dim,
int nstate>
16 template<
int dim,
int nstate>
19 const unsigned int n_dofs_cell = dg->fe_collection[poly_degree].dofs_per_cell;
20 const unsigned int n_quad_pts = dg->volume_quadrature_collection[poly_degree].size();
21 const unsigned int n_shape_fns = n_dofs_cell / nstate;
24 vol_projection.
build_1D_volume_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
27 soln_basis.
build_1D_volume_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
29 dealii::LinearAlgebra::distributed::Vector<double> entropy_var_hat_global(dg->right_hand_side);
30 dealii::LinearAlgebra::distributed::Vector<double> energy_var_hat_global(dg->right_hand_side);
31 std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs_cell);
35 for (
auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
36 if (!cell->is_locally_owned())
continue;
37 cell->get_dof_indices (dofs_indices);
40 std::array<std::vector<double>,nstate> soln_coeff;
41 for(
unsigned int idof=0; idof<n_dofs_cell; idof++){
42 const unsigned int istate = dg->fe_collection[poly_degree].system_to_component_index(idof).first;
43 const unsigned int ishape = dg->fe_collection[poly_degree].system_to_component_index(idof).second;
45 soln_coeff[istate].resize(n_shape_fns);
46 soln_coeff[istate][ishape] = dg->solution(dofs_indices[idof]);
50 std::array<std::vector<double>,nstate> soln_at_q;
51 for(
int istate=0; istate<nstate; istate++){
52 soln_at_q[istate].resize(n_quad_pts);
57 std::array<std::vector<double>,nstate> entropy_var_at_q;
58 std::array<std::vector<double>,nstate> energy_var_at_q;
59 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
60 std::array<double,nstate> soln_state;
61 for(
int istate=0; istate<nstate; istate++){
62 soln_state[istate] = soln_at_q[istate][iquad];
64 std::array<double,nstate> entropy_var_state = euler_double->compute_entropy_variables(soln_state);
65 std::array<double,nstate> kin_energy_state = euler_double->compute_kinetic_energy_variables(soln_state);
66 for(
int istate=0; istate<nstate; istate++){
68 entropy_var_at_q[istate].resize(n_quad_pts);
69 energy_var_at_q[istate].resize(n_quad_pts);
71 energy_var_at_q[istate][iquad] = kin_energy_state[istate];
72 entropy_var_at_q[istate][iquad] = entropy_var_state[istate];
77 for(
int istate=0; istate<nstate; istate++){
79 std::vector<double> entropy_var_hat(n_shape_fns);
82 std::vector<double> energy_var_hat(n_shape_fns);
86 for(
unsigned int ishape=0; ishape<n_shape_fns; ishape++){
87 const unsigned int idof = istate * n_shape_fns + ishape;
88 entropy_var_hat_global[dofs_indices[idof]] = entropy_var_hat[ishape];
89 energy_var_hat_global[dofs_indices[idof]] = energy_var_hat[ishape];
95 dg->assemble_residual();
96 std::array<double,2> change_entropy_and_energy;
97 change_entropy_and_energy[0] = entropy_var_hat_global * dg->right_hand_side;
98 change_entropy_and_energy[1] = energy_var_hat_global * dg->right_hand_side;
99 return change_entropy_and_energy;
101 template<
int dim,
int nstate>
104 const unsigned int n_dofs_cell = dg->fe_collection[poly_degree].dofs_per_cell;
105 const unsigned int n_quad_pts = dg->volume_quadrature_collection[poly_degree].size();
106 const unsigned int n_shape_fns = n_dofs_cell / nstate;
107 const unsigned int grid_degree = dg->high_order_grid->fe_system.tensor_degree();
110 soln_basis.
build_1D_volume_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
111 soln_basis.
build_1D_gradient_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
115 flux_basis.
build_1D_volume_operator(dg->oneD_fe_collection_flux[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
121 flux_basis_stiffness.
build_1D_volume_operator(dg->oneD_fe_collection_flux[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
126 const std::vector<double> &oneD_vol_quad_weights = dg->oneD_quadrature_collection[poly_degree].get_weights();
129 vol_projection.
build_1D_volume_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
131 std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs_cell);
135 double volume_term = 0.0;
136 auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
137 for (
auto cell = dg->dof_handler.begin_active(); cell!= dg->dof_handler.end(); ++cell, ++metric_cell) {
138 if (!cell->is_locally_owned())
continue;
139 cell->get_dof_indices (dofs_indices);
142 const dealii::FESystem<dim> &fe_metric = dg->high_order_grid->fe_system;
143 const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
144 const unsigned int n_grid_nodes = n_metric_dofs / dim;
145 std::vector<dealii::types::global_dof_index> metric_dof_indices(n_metric_dofs);
146 metric_cell->get_dof_indices (metric_dof_indices);
147 std::array<std::vector<double>,dim> mapping_support_points;
148 for(
int idim=0; idim<dim; idim++){
149 mapping_support_points[idim].resize(n_grid_nodes);
153 const std::vector<unsigned int > &index_renumbering = dealii::FETools::hierarchic_to_lexicographic_numbering<dim>(grid_degree);
154 for (
unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
155 const double val = (dg->high_order_grid->volume_nodes[metric_dof_indices[idof]]);
156 const unsigned int istate = fe_metric.system_to_component_index(idof).first;
157 const unsigned int ishape = fe_metric.system_to_component_index(idof).second;
158 const unsigned int igrid_node = index_renumbering[ishape];
159 mapping_support_points[istate][igrid_node] = val;
167 n_quad_pts, n_grid_nodes,
168 mapping_support_points,
170 dg->all_parameters->use_invariant_curl_form);
174 std::array<std::vector<double>,nstate> soln_coeff;
175 for(
unsigned int idof=0; idof<n_dofs_cell; idof++){
176 const unsigned int istate = dg->fe_collection[poly_degree].system_to_component_index(idof).first;
177 const unsigned int ishape = dg->fe_collection[poly_degree].system_to_component_index(idof).second;
179 soln_coeff[istate].resize(n_shape_fns);
180 soln_coeff[istate][ishape] = dg->solution(dofs_indices[idof]);
183 std::array<std::vector<double>,nstate> soln_at_q;
184 std::array<std::vector<double>,dim> vel_at_q;
185 for(
int istate=0; istate<nstate; istate++){
186 soln_at_q[istate].resize(n_quad_pts);
192 std::array<std::vector<double>,nstate> energy_var_vol_int;
193 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
194 std::array<double,nstate> soln_state;
195 for(
int istate=0; istate<nstate; istate++){
196 soln_state[istate] = soln_at_q[istate][iquad];
198 std::array<double,nstate> energy_var;
199 energy_var = euler_double->compute_kinetic_energy_variables(soln_state);
200 for(
int istate=0; istate<nstate; istate++){
202 energy_var_vol_int[istate].resize(n_quad_pts);
204 energy_var_vol_int[istate][iquad] = energy_var[istate];
208 std::array<std::vector<double>,nstate> energy_var_hat;
209 for(
int istate=0; istate<nstate; istate++){
211 energy_var_hat[istate].resize(n_shape_fns);
216 std::array<dealii::Tensor<1,dim,dealii::FullMatrix<double>>,nstate> conv_ref_2pt_flux_at_q;
218 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
220 std::array<double,nstate> soln_state;
221 for(
int istate=0; istate<nstate; istate++){
222 soln_state[istate] = soln_at_q[istate][iquad];
228 dealii::Tensor<2,dim,double> metric_cofactor;
229 for(
int idim=0; idim<dim; idim++){
230 for(
int jdim=0; jdim<dim; jdim++){
238 std::array<dealii::Tensor<1,dim,double>,nstate> conv_phys_flux_2pt;
239 std::vector<std::array<dealii::Tensor<1,dim,double>,nstate>> conv_ref_flux_2pt(n_quad_pts);
240 for (
unsigned int flux_basis=iquad; flux_basis<n_quad_pts; ++flux_basis) {
245 dealii::Tensor<2,dim,double> metric_cofactor_flux_basis;
246 for(
int idim=0; idim<dim; idim++){
247 for(
int jdim=0; jdim<dim; jdim++){
248 metric_cofactor_flux_basis[idim][jdim] = metric_oper.
metric_cofactor_vol[idim][jdim][flux_basis];
251 std::array<double,nstate> soln_state_flux_basis;
252 for(
int istate=0; istate<nstate; istate++){
253 soln_state_flux_basis[istate] = soln_at_q[istate][flux_basis];
256 conv_phys_flux_2pt = euler_double->convective_numerical_split_flux(soln_state, soln_state_flux_basis);
259 const double pressure_int = euler_double->compute_pressure(soln_state);
260 const double pressure_ext = euler_double->compute_pressure(soln_state_flux_basis);
261 for(
int idim=0; idim<dim; idim++){
262 conv_phys_flux_2pt[1+idim][idim] -= 0.5*(pressure_int + pressure_ext);
266 for(
int istate=0; istate<nstate; istate++){
269 conv_phys_flux_2pt[istate],
270 0.5*(metric_cofactor + metric_cofactor_flux_basis),
271 conv_ref_flux_2pt[flux_basis][istate]);
276 for(
int istate=0; istate<nstate; istate++){
280 for(
int idim=0; idim<dim; idim++){
283 conv_ref_2pt_flux_at_q[istate][idim].reinit(n_quad_pts, n_quad_pts);
286 for (
unsigned int flux_basis=iquad; flux_basis<n_quad_pts; ++flux_basis) {
288 conv_ref_2pt_flux_at_q[istate][idim][iquad][flux_basis] = conv_ref_flux_2pt[flux_basis][istate][idim];
289 conv_ref_2pt_flux_at_q[istate][idim][flux_basis][iquad] = conv_ref_flux_2pt[flux_basis][istate][idim];
297 for(
int istate=0; istate<nstate; istate++){
300 std::vector<double> conv_flux_divergence(n_quad_pts);
316 std::vector<double> rhs(n_shape_fns);
319 std::vector<double> ones(n_quad_pts, 1.0);
322 for(
unsigned int ishape=0; ishape<n_shape_fns; ishape++){
323 volume_term += energy_var_hat[istate][ishape] * rhs[ishape];
331 template<
int dim,
int nstate>
334 const unsigned int n_dofs_cell = dg->fe_collection[poly_degree].dofs_per_cell;
335 const unsigned int n_quad_pts = dg->volume_quadrature_collection[poly_degree].size();
336 const unsigned int n_shape_fns = n_dofs_cell / nstate;
340 soln_basis.
build_1D_volume_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
346 std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs_cell);
350 double entropy_fn = 0.0;
351 const std::vector<double> &quad_weights = dg->volume_quadrature_collection[poly_degree].get_weights();
353 auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
354 for (
auto cell = dg->dof_handler.begin_active(); cell!= dg->dof_handler.end(); ++cell, ++metric_cell) {
355 if (!cell->is_locally_owned())
continue;
356 cell->get_dof_indices (dofs_indices);
359 const dealii::FESystem<dim> &fe_metric = dg->high_order_grid->fe_system;
360 const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
361 const unsigned int n_grid_nodes = n_metric_dofs / dim;
362 std::vector<dealii::types::global_dof_index> metric_dof_indices(n_metric_dofs);
363 metric_cell->get_dof_indices (metric_dof_indices);
364 std::array<std::vector<double>,dim> mapping_support_points;
365 for(
int idim=0; idim<dim; idim++){
366 mapping_support_points[idim].resize(n_grid_nodes);
370 const std::vector<unsigned int > &index_renumbering = dealii::FETools::hierarchic_to_lexicographic_numbering<dim>(dg->max_grid_degree);
371 for (
unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
372 const double val = (dg->high_order_grid->volume_nodes[metric_dof_indices[idof]]);
373 const unsigned int istate = fe_metric.system_to_component_index(idof).first;
374 const unsigned int ishape = fe_metric.system_to_component_index(idof).second;
375 const unsigned int igrid_node = index_renumbering[ishape];
376 mapping_support_points[istate][igrid_node] = val;
384 n_quad_pts, n_grid_nodes,
385 mapping_support_points,
387 dg->all_parameters->use_invariant_curl_form);
389 std::array<std::vector<double>,nstate> soln_coeff;
390 for(
unsigned int idof=0; idof<n_dofs_cell; idof++){
391 const unsigned int istate = dg->fe_collection[poly_degree].system_to_component_index(idof).first;
392 const unsigned int ishape = dg->fe_collection[poly_degree].system_to_component_index(idof).second;
394 soln_coeff[istate].resize(n_shape_fns);
395 soln_coeff[istate][ishape] = dg->solution(dofs_indices[idof]);
398 std::array<std::vector<double>,nstate> soln_at_q;
399 for(
int istate=0; istate<nstate; istate++){
400 soln_at_q[istate].resize(n_quad_pts);
404 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
405 std::array<double,nstate> soln_state;
406 for(
int istate=0; istate<nstate; istate++){
407 soln_state[istate] = soln_at_q[istate][iquad];
409 const double density = soln_state[0];
410 const double pressure = euler_double->compute_pressure(soln_state);
411 const double entropy = log(pressure) - euler_double->gam * log(density);
412 const double quadrature_entropy = -density*entropy/euler_double->gamm1;
414 entropy_fn += quadrature_entropy * quad_weights[iquad] * metric_oper.
det_Jac_vol[iquad];
422 template<
int dim,
int nstate>
426 int overintegrate = 10 ;
427 dealii::QGauss<1> quad_extra(dg->max_degree+1+overintegrate);
428 const unsigned int n_dofs_cell = dg->fe_collection[poly_degree].dofs_per_cell;
429 const unsigned int n_quad_pts = dg->volume_quadrature_collection[poly_degree].size();
430 const unsigned int n_shape_fns = n_dofs_cell / nstate;
434 soln_basis.
build_1D_volume_operator(dg->oneD_fe_collection_1state[poly_degree], dg->oneD_quadrature_collection[poly_degree]);
440 std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs_cell);
442 double total_kinetic_energy = 0;
444 const std::vector<double> &quad_weights = dg->volume_quadrature_collection[poly_degree].get_weights();
446 auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
447 for (
auto cell = dg->dof_handler.begin_active(); cell!= dg->dof_handler.end(); ++cell, ++metric_cell) {
448 if (!cell->is_locally_owned())
continue;
451 cell->get_dof_indices (dofs_indices);
453 const dealii::FESystem<dim> &fe_metric = dg->high_order_grid->fe_system;
454 const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
455 const unsigned int n_grid_nodes = n_metric_dofs / dim;
456 std::vector<dealii::types::global_dof_index> metric_dof_indices(n_metric_dofs);
457 metric_cell->get_dof_indices (metric_dof_indices);
458 std::array<std::vector<double>,dim> mapping_support_points;
459 for(
int idim=0; idim<dim; idim++){
460 mapping_support_points[idim].resize(n_grid_nodes);
464 const std::vector<unsigned int > &index_renumbering = dealii::FETools::hierarchic_to_lexicographic_numbering<dim>(dg->max_grid_degree);
465 for (
unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
466 const double val = (dg->high_order_grid->volume_nodes[metric_dof_indices[idof]]);
467 const unsigned int istate = fe_metric.system_to_component_index(idof).first;
468 const unsigned int ishape = fe_metric.system_to_component_index(idof).second;
469 const unsigned int igrid_node = index_renumbering[ishape];
470 mapping_support_points[istate][igrid_node] = val;
478 n_quad_pts, n_grid_nodes,
479 mapping_support_points,
481 dg->all_parameters->use_invariant_curl_form);
483 std::array<std::vector<double>,nstate> soln_coeff;
484 for(
unsigned int idof=0; idof<n_dofs_cell; idof++){
485 const unsigned int istate = dg->fe_collection[poly_degree].system_to_component_index(idof).first;
486 const unsigned int ishape = dg->fe_collection[poly_degree].system_to_component_index(idof).second;
488 soln_coeff[istate].resize(n_shape_fns);
489 soln_coeff[istate][ishape] = dg->solution(dofs_indices[idof]);
492 std::array<std::vector<double>,nstate> soln_at_q;
493 for(
int istate=0; istate<nstate; istate++){
494 soln_at_q[istate].resize(n_quad_pts);
498 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
499 std::array<double,nstate> soln_state;
500 for(
int istate=0; istate<nstate; istate++){
501 soln_state[istate] = soln_at_q[istate][iquad];
503 const double density = soln_state[0];
504 const double quadrature_kinetic_energy = 0.5*(soln_state[1]*soln_state[1]
505 + soln_state[2]*soln_state[2]
506 + soln_state[3]*soln_state[3])/density;
508 total_kinetic_energy += quadrature_kinetic_energy * quad_weights[iquad] * metric_oper.
det_Jac_vol[iquad];
511 return total_kinetic_energy;
514 template<
int dim,
int nstate>
518 const unsigned int n_dofs_cell = nstate*pow(poly_degree+1,dim);
519 const unsigned int n_quad_pts = pow(poly_degree+1,dim);
520 std::vector<dealii::types::global_dof_index> dofs_indices1 (n_dofs_cell);
522 double cfl_min = 1e100;
524 for (
auto cell = dg->dof_handler.begin_active(); cell!=dg->dof_handler.end(); ++cell) {
525 if (!cell->is_locally_owned())
continue;
527 cell->get_dof_indices (dofs_indices1);
528 std::vector< std::array<double,nstate>> soln_at_q(n_quad_pts);
529 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
530 for (
int istate=0; istate<nstate; istate++) {
531 soln_at_q[iquad][istate] = 0;
534 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
535 dealii::Point<dim> qpoint = dg->volume_quadrature_collection[poly_degree].point(iquad);
536 for(
unsigned int idof=0; idof<n_dofs_cell; idof++){
537 const unsigned int istate = dg->fe_collection[poly_degree].system_to_component_index(idof).first;
538 soln_at_q[iquad][istate] += dg->solution[dofs_indices1[idof]] * dg->fe_collection[poly_degree].shape_value_component(idof, qpoint, istate);
542 std::vector< double > convective_eigenvalues(n_quad_pts);
543 for (
unsigned int isol = 0; isol < n_quad_pts; ++isol) {
544 convective_eigenvalues[isol] = pde_physics_double->max_convective_eigenvalue (soln_at_q[isol]);
546 const double max_eig = *(std::max_element(convective_eigenvalues.begin(), convective_eigenvalues.end()));
548 double cfl = 0.1 * delta_x/max_eig;
556 template <
int dim,
int nstate>
559 using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
560 std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
562 typename dealii::Triangulation<dim>::MeshSmoothing(
563 dealii::Triangulation<dim>::smoothing_on_refinement |
564 dealii::Triangulation<dim>::smoothing_on_coarsening));
570 double right = 2 * dealii::numbers::PI;
571 const int n_refinements = 2;
572 unsigned int poly_degree = 3;
577 PHiLiP::Grids::nonsymmetric_curved_grid<dim,Triangulation>(*grid, n_refinements);
581 PHiLiP::Grids::straight_periodic_cube<dim,Triangulation>(grid, left, right, pow(2.0,n_refinements));
586 dg->allocate_system (
false,
false,
false);
588 pcout <<
"Implement initial conditions" << std::endl;
589 std::shared_ptr< InitialConditionFunction<dim,nstate,double> > initial_condition_function =
593 const unsigned int n_global_active_cells2 = grid->n_global_active_cells();
594 double delta_x = (right-left)/pow(n_global_active_cells2,1.0/dim)/(poly_degree+1.0);
595 pcout<<
" delta x "<<delta_x<<std::endl;
599 pcout <<
"creating ODE solver" << std::endl;
601 pcout <<
"ODE solver successfully created" << std::endl;
602 double finalTime = 14.;
604 pcout <<
" number dofs " << dg->dof_handler.n_dofs()<<std::endl;
605 pcout <<
"preparing to advance solution in time" << std::endl;
607 ode_solver->current_iteration = 0;
608 ode_solver->allocate_ode_system();
611 const double initial_energy_mpi = (dealii::Utilities::MPI::sum(initial_energy,
mpi_communicator));
613 const double initial_entropy_mpi = (dealii::Utilities::MPI::sum(initial_entropy,
mpi_communicator));
615 std::ofstream myfile (all_parameters_new.
energy_file +
".gpl" , std::ios::trunc);
617 while(ode_solver->current_time < finalTime){
619 const double time_step =
get_timestep(dg,poly_degree, delta_x);
621 pcout<<
"time step "<<time_step<<
" current time "<<ode_solver->current_time<<std::endl;
625 ode_solver->step_in_time(dt,
false);
626 ode_solver->current_iteration += 1;
629 if (is_output_iteration) {
631 dg->output_results_vtk(file_number);
636 const double current_change_entropy_mpi = dealii::Utilities::MPI::sum(current_change_entropy[0],
mpi_communicator);
637 const double current_change_energy_mpi = dealii::Utilities::MPI::sum(current_change_entropy[1],
mpi_communicator);
639 myfile<<ode_solver->current_time<<
" "<< current_change_entropy_mpi <<std::endl;
640 pcout <<
"M plus K norm Change in Entropy at time " << ode_solver->current_time <<
" is " << current_change_entropy_mpi<< std::endl;
641 pcout <<
"M plus K norm Change in Kinetic Energy at time " << ode_solver->current_time <<
" is " << current_change_energy_mpi<< std::endl;
643 if(abs(current_change_entropy[0]) > 1e-12 && (dg->all_parameters->two_point_num_flux_type == Parameters::AllParameters::TwoPointNumericalFlux::IR || dg->all_parameters->two_point_num_flux_type == Parameters::AllParameters::TwoPointNumericalFlux::CH || dg->all_parameters->two_point_num_flux_type == Parameters::AllParameters::TwoPointNumericalFlux::Ra)){
644 pcout <<
" Change in entropy was not monotonically conserved." << std::endl;
650 const double current_energy_mpi = (dealii::Utilities::MPI::sum(current_energy,
mpi_communicator));
651 pcout <<
"Normalized kinetic energy " << ode_solver->current_time <<
" is " << current_energy_mpi/initial_energy_mpi<< std::endl;
654 const double current_entropy_mpi = (dealii::Utilities::MPI::sum(current_entropy,
mpi_communicator));
655 pcout <<
"Normalized entropy " << ode_solver->current_time <<
" is " << current_entropy_mpi/initial_entropy_mpi<< std::endl;
659 double current_vol_work_mpi = (dealii::Utilities::MPI::sum(current_vol_work,
mpi_communicator));
660 pcout<<
"volume work "<<current_vol_work_mpi<<std::endl;
661 myfile << ode_solver->current_time <<
" " << std::fixed << std::setprecision(16) << current_vol_work_mpi<< std::endl;
664 if(abs(current_vol_work_mpi) > 1e-12 && (dg->all_parameters->two_point_num_flux_type == Parameters::AllParameters::TwoPointNumericalFlux::Ra || dg->all_parameters->two_point_num_flux_type == Parameters::AllParameters::TwoPointNumericalFlux::KG ) ){
665 pcout<<
"The kinetic energy volume work is not zero."<<std::endl;
dealii::Tensor< 2, dim, std::vector< real > > metric_cofactor_vol
The volume metric cofactor matrix.
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
bool use_curvilinear_grid
Flag to use curvilinear grid.
void divergence_two_pt_flux_Hadamard_product(const dealii::Tensor< 1, dim, dealii::FullMatrix< real >> &input_mat, std::vector< real > &output_vect, const std::vector< real > &weights, const dealii::FullMatrix< double > &basis, const double scaling=2.0)
Computes the divergence of the 2pt flux Hadamard products, then sums the rows.
double compute_kinetic_energy(const std::shared_ptr< DGBase< dim, double > > &dg, unsigned int poly_degree) const
Computes kinetic energy.
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.
int output_solution_every_x_steps
Outputs the solution every x steps to .vtk file.
void inner_product_1D(const std::vector< real > &input_vect, const std::vector< real > &weight_vect, std::vector< real > &output_vect, const dealii::FullMatrix< double > &basis_x, const bool adding=false, const double factor=1.0)
Apply the inner product operation using the 1D operator in each direction.
Files for the baseline physics.
unsigned int print_iteration_modulo
If ode_output==verbose, print every print_iteration_modulo iterations.
dealii::FullMatrix< double > oneD_skew_symm_vol_oper
Skew-symmetric volume operator .
void build_1D_shape_functions_at_grid_nodes(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Constructs the volume operator and gradient operator.
double compute_volume_term(const std::shared_ptr< DGBase< dim, double > > &dg, unsigned int poly_degree) const
Computes the volume term kinetic energy production.
Main parameter class that contains the various other sub-parameter classes.
std::string energy_file
Energy file.
std::vector< real > det_Jac_vol
The determinant of the metric Jacobian at volume cubature nodes.
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.
void build_1D_gradient_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
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.
double get_timestep(const std::shared_ptr< DGBase< dim, double > > &dg, unsigned int poly_degree, const double delta_x) const
Computes the timestep from max eignevector.
Base metric operators class that stores functions used in both the volume and on surface.
The mapping shape functions evaluated at the desired nodes (facet set included in volume grid nodes f...
EulerTaylorGreen(const Parameters::AllParameters *const parameters_input)
Constructor.
Euler equations. Derived from PhysicsBase.
Euler Taylor Green Vortex.
void transform_physical_to_reference(const dealii::Tensor< 1, dim, real > &phys, const dealii::Tensor< 2, dim, real > &metric_cofactor, dealii::Tensor< 1, dim, real > &ref)
Given a physical tensor, return the reference tensor.
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
static std::shared_ptr< PhysicsBase< dim, nstate, real > > create_Physics(const Parameters::AllParameters *const parameters_input, std::shared_ptr< ModelBase< dim, nstate, real > > model_input=nullptr)
Factory to return the correct physics given input file.
void build_1D_surface_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 0 > &quadrature)
Assembles the one dimensional operator.
void build_volume_metric_operators(const unsigned int n_quad_pts, const unsigned int n_metric_dofs, const std::array< std::vector< real >, dim > &mapping_support_points, mapping_shape_functions< dim, n_faces, real > &mapping_basis, const bool use_invariant_curl_form=false)
Builds the volume metric operators.
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
dealii::ConditionalOStream pcout
ConditionalOStream.
DGBase is independent of the number of state variables.
double compute_entropy(const std::shared_ptr< DGBase< dim, double > > &dg, unsigned int poly_degree) const
Computes entropy in the norm.
void build_1D_shape_functions_at_flux_nodes(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature, const dealii::Quadrature< 0 > &face_quadrature)
Constructs the volume, gradient, surface, and surface gradient operator.
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.
Local stiffness matrix without jacobian dependence.
Base class of all the tests.
int overintegration
Number of additional quadrature points to use.
int run_test() const override
Ensure that the kinetic energy is bounded.
static std::shared_ptr< InitialConditionFunction< dim, nstate, real > > create_InitialConditionFunction(Parameters::AllParameters const *const param)
Construct InitialConditionFunction object from global parameter file.
std::array< double, 2 > compute_change_in_entropy(const std::shared_ptr< DGBase< dim, double > > &dg, unsigned int poly_degree) const
Computes change in entropy in the norm.