25#ifndef EXADG_OPERATORS_ADAPTIVE_MESH_REFINEMENT_H_
26#define EXADG_OPERATORS_ADAPTIVE_MESH_REFINEMENT_H_
29#include <deal.II/base/quadrature_lib.h>
30#include <deal.II/distributed/grid_refinement.h>
31#include <deal.II/dofs/dof_tools.h>
32#include <deal.II/grid/grid_refinement.h>
33#include <deal.II/numerics/error_estimator.h>
34#include <deal.II/numerics/solution_transfer.h>
37#include <exadg/utilities/print_functions.h>
41struct AdaptiveMeshRefinementData
43 AdaptiveMeshRefinementData()
44 : trigger_every_n_time_steps(1),
45 maximum_refinement_level(10),
46 minimum_refinement_level(0),
47 preserve_boundary_cells(
false),
48 fraction_of_cells_to_be_refined(0.0),
49 fraction_of_cells_to_be_coarsened(0.0)
54 print(dealii::ConditionalOStream
const & pcout)
const
56 print_parameter(pcout,
"Enable adaptivity",
true);
57 print_parameter(pcout,
"Triggered every n time steps", trigger_every_n_time_steps);
58 print_parameter(pcout,
"Maximum refinement level", maximum_refinement_level);
59 print_parameter(pcout,
"Minimum refinement level", minimum_refinement_level);
60 print_parameter(pcout,
"Preserve boundary cells", preserve_boundary_cells);
61 print_parameter(pcout,
"Fraction of cells to be refined", fraction_of_cells_to_be_refined);
62 print_parameter(pcout,
"Fraction of cells to be coarsened", fraction_of_cells_to_be_coarsened);
65 unsigned int trigger_every_n_time_steps;
67 int maximum_refinement_level;
68 int minimum_refinement_level;
69 bool preserve_boundary_cells;
71 double fraction_of_cells_to_be_refined;
72 double fraction_of_cells_to_be_coarsened;
80 unsigned int const time_step_number)
82 return ((time_step_number == 0) or not(time_step_number % trigger_every_n_time_steps));
95 for(
auto & cell : tria.active_cell_iterators())
97 if(cell->is_locally_owned())
100 if(cell->level() == amr_data.maximum_refinement_level)
102 cell->clear_refine_flag();
106 if(cell->level() <= amr_data.minimum_refinement_level)
108 cell->clear_coarsen_flag();
112 if(amr_data.preserve_boundary_cells and cell->at_boundary())
114 cell->clear_refine_flag();
115 cell->clear_coarsen_flag();
128 bool any_flag_set =
false;
130 for(
auto & cell : tria.active_cell_iterators())
132 if(cell->is_locally_owned())
134 if(cell->refine_flag_set() or cell->coarsen_flag_set())
142 any_flag_set = dealii::Utilities::MPI::logical_or(any_flag_set, tria.get_mpi_communicator());
147template<
int dim,
typename Number,
typename VectorType>
149mark_cells_kelly_error_estimator(dealii::Triangulation<dim> & tria,
150 dealii::DoFHandler<dim>
const & dof_handler,
151 dealii::AffineConstraints<Number>
const & constraints,
152 dealii::Mapping<dim>
const & mapping,
153 VectorType
const & solution,
154 unsigned int const n_face_quadrature_points,
155 AdaptiveMeshRefinementData
const & amr_data)
157 VectorType locally_relevant_solution;
158 locally_relevant_solution.reinit(dof_handler.locally_owned_dofs(),
159 dealii::DoFTools::extract_locally_relevant_dofs(dof_handler),
160 dof_handler.get_mpi_communicator());
161 locally_relevant_solution.copy_locally_owned_data_from(solution);
162 constraints.distribute(locally_relevant_solution);
163 locally_relevant_solution.update_ghost_values();
165 dealii::QGauss<dim - 1> face_quadrature(n_face_quadrature_points);
167 dealii::Vector<float> estimated_error_per_cell(tria.n_active_cells());
169 dealii::KellyErrorEstimator<dim>::estimate(mapping,
173 locally_relevant_solution,
174 estimated_error_per_cell);
176 dealii::parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number(
178 estimated_error_per_cell,
179 amr_data.fraction_of_cells_to_be_refined,
180 amr_data.fraction_of_cells_to_be_coarsened);
void limit_coarsening_and_refinement(dealii::Triangulation< dim > &tria, AdaptiveMeshRefinementData const &amr_data)
Definition adaptive_mesh_refinement.h:92
bool any_cells_flagged_for_coarsening_or_refinement(dealii::Triangulation< dim > const &tria)
Definition adaptive_mesh_refinement.h:126
bool trigger_coarsening_and_refinement_now(unsigned int const trigger_every_n_time_steps, unsigned int const time_step_number)
Definition adaptive_mesh_refinement.h:79
Definition adaptive_mesh_refinement.h:42