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/distributed/solution_transfer.h>
32#include <deal.II/dofs/dof_tools.h>
33#include <deal.II/grid/grid_refinement.h>
34#include <deal.II/numerics/error_estimator.h>
35#include <deal.II/numerics/solution_transfer.h>
38#include <exadg/utilities/print_functions.h>
42struct AdaptiveMeshRefinementData
44 AdaptiveMeshRefinementData()
45 : trigger_every_n_time_steps(1),
46 maximum_refinement_level(10),
47 minimum_refinement_level(0),
48 preserve_boundary_cells(
false),
49 fraction_of_cells_to_be_refined(0.0),
50 fraction_of_cells_to_be_coarsened(0.0)
55 print(dealii::ConditionalOStream
const & pcout)
const
57 print_parameter(pcout,
"Enable adaptivity",
true);
58 print_parameter(pcout,
"Triggered every n time steps", trigger_every_n_time_steps);
59 print_parameter(pcout,
"Maximum refinement level", maximum_refinement_level);
60 print_parameter(pcout,
"Minimum refinement level", minimum_refinement_level);
61 print_parameter(pcout,
"Preserve boundary cells", preserve_boundary_cells);
62 print_parameter(pcout,
"Fraction of cells to be refined", fraction_of_cells_to_be_refined);
63 print_parameter(pcout,
"Fraction of cells to be coarsened", fraction_of_cells_to_be_coarsened);
66 unsigned int trigger_every_n_time_steps;
68 int maximum_refinement_level;
69 int minimum_refinement_level;
70 bool preserve_boundary_cells;
72 double fraction_of_cells_to_be_refined;
73 double fraction_of_cells_to_be_coarsened;
81 unsigned int const time_step_number)
83 return ((time_step_number == 0) or not(time_step_number % trigger_every_n_time_steps));
96 for(
auto & cell : tria.active_cell_iterators())
98 if(cell->is_locally_owned())
101 if(cell->level() == amr_data.maximum_refinement_level)
103 cell->clear_refine_flag();
107 if(cell->level() <= amr_data.minimum_refinement_level)
109 cell->clear_coarsen_flag();
113 if(amr_data.preserve_boundary_cells and cell->at_boundary())
115 cell->clear_refine_flag();
116 cell->clear_coarsen_flag();
129 bool any_flag_set =
false;
131 for(
auto & cell : tria.active_cell_iterators())
133 if(cell->is_locally_owned())
135 if(cell->refine_flag_set() or cell->coarsen_flag_set())
143 any_flag_set = dealii::Utilities::MPI::logical_or(any_flag_set, tria.get_communicator());
148template<
int dim,
typename Number,
typename VectorType>
150mark_cells_kelly_error_estimator(dealii::Triangulation<dim> & tria,
151 dealii::DoFHandler<dim>
const & dof_handler,
152 dealii::AffineConstraints<Number>
const & constraints,
153 dealii::Mapping<dim>
const & mapping,
154 VectorType
const & solution,
155 unsigned int const n_face_quadrature_points,
156 AdaptiveMeshRefinementData
const & amr_data)
158 VectorType locally_relevant_solution;
159 locally_relevant_solution.reinit(dof_handler.locally_owned_dofs(),
160 dealii::DoFTools::extract_locally_relevant_dofs(dof_handler),
161 dof_handler.get_communicator());
162 locally_relevant_solution.copy_locally_owned_data_from(solution);
163 constraints.distribute(locally_relevant_solution);
164 locally_relevant_solution.update_ghost_values();
166 dealii::QGauss<dim - 1> face_quadrature(n_face_quadrature_points);
168 dealii::Vector<float> estimated_error_per_cell(tria.n_active_cells());
170 dealii::KellyErrorEstimator<dim>::estimate(mapping,
174 locally_relevant_solution,
175 estimated_error_per_cell);
177 dealii::parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number(
179 estimated_error_per_cell,
180 amr_data.fraction_of_cells_to_be_refined,
181 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:93
bool any_cells_flagged_for_coarsening_or_refinement(dealii::Triangulation< dim > const &tria)
Definition adaptive_mesh_refinement.h:127
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:80
Definition adaptive_mesh_refinement.h:43