[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
grid_refinement_uniform.cpp
1 #include "grid_refinement_uniform.h"
2 
3 #include <deal.II/dofs/dof_tools.h>
4 
5 #include "dg/dg_base.hpp"
6 #include "mesh/high_order_grid.h"
7 #include "parameters/parameters_grid_refinement.h"
8 
9 namespace PHiLiP {
10 
11 namespace GridRefinement {
12 
13 template <int dim, int nstate, typename real, typename MeshType>
15 {
17  RefinementTypeEnum refinement_type = this->grid_refinement_param.refinement_type;
18 
19  // setting up the solution transfer
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);
23 
24  dealii::LinearAlgebra::distributed::Vector<real> solution_old(this->dg->solution);
25  solution_old.update_ghost_values();
26 
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);
31 
32  this->dg->high_order_grid->prepare_for_coarsening_and_refinement();
33  this->dg->triangulation->prepare_coarsening_and_refinement();
34 
35  if(refinement_type == RefinementTypeEnum::h){
36  refine_grid_h();
37  }else if(refinement_type == RefinementTypeEnum::p){
38  refine_grid_p();
39  }else if(refinement_type == RefinementTypeEnum::hp){
40  refine_grid_hp();
41  }
42 
43  this->tria->execute_coarsening_and_refinement();
44  this->dg->high_order_grid->execute_coarsening_and_refinement();
45 
46  // transfering the solution from solution_old
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();
51 
52  // increase the count
53  this->iteration++;
54 }
55 
56 // functions for the refinement calls for each of the classes
57 template <int dim, int nstate, typename real, typename MeshType>
59 {
60  this->tria->set_all_refine_flags();
61 }
62 
63 template <int dim, int nstate, typename real, typename MeshType>
65 {
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);
69 
70 }
71 
72 template <int dim, int nstate, typename real, typename MeshType>
74 {
75  refine_grid_h();
76  refine_grid_p();
77 }
78 
79 template <int dim, int nstate, typename real, typename MeshType>
80 std::vector< std::pair<dealii::Vector<real>, std::string> > GridRefinement_Uniform<dim,nstate,real,MeshType>::output_results_vtk_method()
81 {
82  // nothing special to do here
83  std::vector< std::pair<dealii::Vector<real>, std::string> > data_out_vector;
84 
85  return data_out_vector;
86 }
87 
88 // dealii::Triangulation<PHILIP_DIM>
94 
95 // dealii::parallel::shared::Triangulation<PHILIP_DIM>
101 
102 #if PHILIP_DIM != 1
103 // dealii::parallel::distributed::Triangulation<PHILIP_DIM>
109 #endif
110 
111 } // namespace GridRefinement
112 
113 } // namespace PHiLiP
Files for the baseline physics.
Definition: ADTypes.hpp:10
RefinementType
Controls the type of refinement to be performed.
void refine_grid() override
Perform call to the grid refinement object of choice.
std::vector< std::pair< dealii::Vector< real >, std::string > > output_results_vtk_method() override
Output refinement method dependent results.