1 #include "grid_refinement_uniform.h" 3 #include <deal.II/dofs/dof_tools.h> 5 #include "dg/dg_base.hpp" 6 #include "mesh/high_order_grid.h" 7 #include "parameters/parameters_grid_refinement.h" 11 namespace GridRefinement {
13 template <
int dim,
int nstate,
typename real,
typename MeshType>
17 RefinementTypeEnum refinement_type = this->grid_refinement_param.refinement_type;
20 dealii::IndexSet locally_owned_dofs, locally_relevant_dofs;
21 locally_owned_dofs = this->dg->dof_handler.locally_owned_dofs();
22 dealii::DoFTools::extract_locally_relevant_dofs(this->dg->dof_handler, locally_relevant_dofs);
24 dealii::LinearAlgebra::distributed::Vector<real> solution_old(this->dg->solution);
25 solution_old.update_ghost_values();
27 dealii::parallel::distributed::SolutionTransfer<
28 dim, dealii::LinearAlgebra::distributed::Vector<real>, dealii::DoFHandler<dim>
29 > solution_transfer(this->dg->dof_handler);
30 solution_transfer.prepare_for_coarsening_and_refinement(solution_old);
32 this->dg->high_order_grid->prepare_for_coarsening_and_refinement();
33 this->dg->triangulation->prepare_coarsening_and_refinement();
35 if(refinement_type == RefinementTypeEnum::h){
37 }
else if(refinement_type == RefinementTypeEnum::p){
39 }
else if(refinement_type == RefinementTypeEnum::hp){
43 this->tria->execute_coarsening_and_refinement();
44 this->dg->high_order_grid->execute_coarsening_and_refinement();
47 this->dg->allocate_system();
48 this->dg->solution.zero_out_ghosts();
49 solution_transfer.interpolate(this->dg->solution);
50 this->dg->solution.update_ghost_values();
57 template <
int dim,
int nstate,
typename real,
typename MeshType>
60 this->tria->set_all_refine_flags();
63 template <
int dim,
int nstate,
typename real,
typename MeshType>
66 for(
auto cell = this->dg->dof_handler.begin_active(); cell != this->dg->dof_handler.end(); ++cell)
67 if(cell->is_locally_owned() && cell->active_fe_index()+1 <= this->dg->max_degree)
68 cell->set_future_fe_index(cell->active_fe_index()+1);
72 template <
int dim,
int nstate,
typename real,
typename MeshType>
79 template <
int dim,
int nstate,
typename real,
typename MeshType>
83 std::vector< std::pair<dealii::Vector<real>, std::string> > data_out_vector;
85 return data_out_vector;
Files for the baseline physics.
RefinementType
Controls the type of refinement to be performed.