ExaDG
Loading...
Searching...
No Matches
adaptive_mesh_refinement.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 *
21 * Source: https://github.com/MeltPoolDG/MeltPoolDG
22 * Author: Peter Munch, Magdalena Schreter, TUM, December 2020
23 */
24
25#ifndef EXADG_OPERATORS_ADAPTIVE_MESH_REFINEMENT_H_
26#define EXADG_OPERATORS_ADAPTIVE_MESH_REFINEMENT_H_
27
28// deal.II
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>
35
36// ExaDG
37#include <exadg/utilities/print_functions.h>
38
39namespace ExaDG
40{
41struct AdaptiveMeshRefinementData
42{
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)
50 {
51 }
52
53 void
54 print(dealii::ConditionalOStream const & pcout) const
55 {
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);
63 }
64
65 unsigned int trigger_every_n_time_steps;
66
67 int maximum_refinement_level;
68 int minimum_refinement_level;
69 bool preserve_boundary_cells;
70
71 double fraction_of_cells_to_be_refined;
72 double fraction_of_cells_to_be_coarsened;
73};
74
78inline bool
79trigger_coarsening_and_refinement_now(unsigned int const trigger_every_n_time_steps,
80 unsigned int const time_step_number)
81{
82 return ((time_step_number == 0) or not(time_step_number % trigger_every_n_time_steps));
83}
84
90template<int dim>
91void
92limit_coarsening_and_refinement(dealii::Triangulation<dim> & tria,
93 AdaptiveMeshRefinementData const & amr_data)
94{
95 for(auto & cell : tria.active_cell_iterators())
96 {
97 if(cell->is_locally_owned())
98 {
99 // Clear refinement flags on maximum refinement level.
100 if(cell->level() == amr_data.maximum_refinement_level)
101 {
102 cell->clear_refine_flag();
103 }
104
105 // Clear coarsening flags on minimum refinement level.
106 if(cell->level() <= amr_data.minimum_refinement_level)
107 {
108 cell->clear_coarsen_flag();
109 }
110
111 // Do not coarsen/refine cells along boundary if requested
112 if(amr_data.preserve_boundary_cells and cell->at_boundary())
113 {
114 cell->clear_refine_flag();
115 cell->clear_coarsen_flag();
116 }
117 }
118 }
119}
120
124template<int dim>
125bool
126any_cells_flagged_for_coarsening_or_refinement(dealii::Triangulation<dim> const & tria)
127{
128 bool any_flag_set = false;
129
130 for(auto & cell : tria.active_cell_iterators())
131 {
132 if(cell->is_locally_owned())
133 {
134 if(cell->refine_flag_set() or cell->coarsen_flag_set())
135 {
136 any_flag_set = true;
137 break;
138 }
139 }
140 }
141
142 any_flag_set = dealii::Utilities::MPI::logical_or(any_flag_set, tria.get_mpi_communicator());
143
144 return any_flag_set;
145}
146
147template<int dim, typename Number, typename VectorType>
148void
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)
156{
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();
164
165 dealii::QGauss<dim - 1> face_quadrature(n_face_quadrature_points);
166
167 dealii::Vector<float> estimated_error_per_cell(tria.n_active_cells());
168
169 dealii::KellyErrorEstimator<dim>::estimate(mapping,
170 dof_handler,
171 face_quadrature,
172 {}, // Neumann BC
173 locally_relevant_solution,
174 estimated_error_per_cell);
175
176 dealii::parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number(
177 tria,
178 estimated_error_per_cell,
179 amr_data.fraction_of_cells_to_be_refined,
180 amr_data.fraction_of_cells_to_be_coarsened);
181}
182
183} // namespace ExaDG
184
185#endif /* EXADG_OPERATORS_ADAPTIVE_MESH_REFINEMENT_H_ */
Definition driver.cpp:33
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