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 INCLUDE_EXADG_OPERATORS_ADAPTIVE_MESH_REFINEMENT_H_
26#define INCLUDE_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/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>
36
37// ExaDG
38#include <exadg/utilities/print_functions.h>
39
40
41namespace ExaDG
42{
44{
46 : trigger_every_n_time_steps(1),
47 maximum_refinement_level(10),
48 minimum_refinement_level(0),
49 preserve_boundary_cells(false),
50 fraction_of_cells_to_be_refined(0.0),
51 fraction_of_cells_to_be_coarsened(0.0)
52 {
53 }
54
55 void
56 print(dealii::ConditionalOStream const & pcout) const
57 {
58 print_parameter(pcout, "Enable adaptivity", true);
59 print_parameter(pcout, "Triggered every n time steps", trigger_every_n_time_steps);
60 print_parameter(pcout, "Maximum refinement level", maximum_refinement_level);
61 print_parameter(pcout, "Minimum refinement level", minimum_refinement_level);
62 print_parameter(pcout, "Preserve boundary cells", preserve_boundary_cells);
63 print_parameter(pcout, "Fraction of cells to be refined", fraction_of_cells_to_be_refined);
64 print_parameter(pcout, "Fraction of cells to be coarsened", fraction_of_cells_to_be_coarsened);
65 }
66
67 unsigned int trigger_every_n_time_steps;
68
69 int maximum_refinement_level;
70 int minimum_refinement_level;
71 bool preserve_boundary_cells;
72
73 double fraction_of_cells_to_be_refined;
74 double fraction_of_cells_to_be_coarsened;
75};
76
80inline bool
81trigger_coarsening_and_refinement_now(unsigned int const trigger_every_n_time_steps,
82 unsigned int const time_step_number)
83{
84 return ((time_step_number == 0) or not(time_step_number % trigger_every_n_time_steps));
85}
86
92template<int dim>
93void
94limit_coarsening_and_refinement(dealii::Triangulation<dim> & tria,
95 AdaptiveMeshRefinementData const & amr_data)
96{
97 for(auto & cell : tria.active_cell_iterators())
98 {
99 if(cell->is_locally_owned())
100 {
101 // Clear refinement flags on maximum refinement level.
102 if(cell->level() == amr_data.maximum_refinement_level)
103 {
104 cell->clear_refine_flag();
105 }
106
107 // Clear coarsening flags on minimum refinement level.
108 if(cell->level() <= amr_data.minimum_refinement_level)
109 {
110 cell->clear_coarsen_flag();
111 }
112
113 // Do not coarsen/refine cells along boundary if requested
114 if(amr_data.preserve_boundary_cells and cell->at_boundary())
115 {
116 cell->clear_refine_flag();
117 cell->clear_coarsen_flag();
118 }
119 }
120 }
121}
122
126template<int dim>
127bool
128any_cells_flagged_for_coarsening_or_refinement(dealii::Triangulation<dim> const & tria)
129{
130 bool any_flag_set = false;
131
132 for(auto & cell : tria.active_cell_iterators())
133 {
134 if(cell->is_locally_owned())
135 {
136 if(cell->refine_flag_set() or cell->coarsen_flag_set())
137 {
138 any_flag_set = true;
139 break;
140 }
141 }
142 }
143
144 any_flag_set = dealii::Utilities::MPI::logical_or(any_flag_set, tria.get_communicator());
145
146 return any_flag_set;
147}
148
149template<int dim, typename Number, typename VectorType>
150void
151mark_cells_kelly_error_estimator(dealii::Triangulation<dim> & tria,
152 dealii::DoFHandler<dim> const & dof_handler,
153 dealii::AffineConstraints<Number> const & constraints,
154 dealii::Mapping<dim> const & mapping,
155 VectorType const & solution,
156 unsigned int const n_face_quadrature_points,
157 AdaptiveMeshRefinementData const & amr_data)
158{
159 VectorType locally_relevant_solution;
160 locally_relevant_solution.reinit(dof_handler.locally_owned_dofs(),
161 dealii::DoFTools::extract_locally_relevant_dofs(dof_handler),
162 dof_handler.get_communicator());
163 locally_relevant_solution.copy_locally_owned_data_from(solution);
164 constraints.distribute(locally_relevant_solution);
165 locally_relevant_solution.update_ghost_values();
166
167 dealii::QGauss<dim - 1> face_quadrature(n_face_quadrature_points);
168
169 dealii::Vector<float> estimated_error_per_cell(tria.n_active_cells());
170
171 dealii::KellyErrorEstimator<dim>::estimate(mapping,
172 dof_handler,
173 face_quadrature,
174 {}, // Neumann BC
175 locally_relevant_solution,
176 estimated_error_per_cell);
177
178 dealii::parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number(
179 tria,
180 estimated_error_per_cell,
181 amr_data.fraction_of_cells_to_be_refined,
182 amr_data.fraction_of_cells_to_be_coarsened);
183}
184
185
186} // namespace ExaDG
187
188#endif /* INCLUDE_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:94
bool any_cells_flagged_for_coarsening_or_refinement(dealii::Triangulation< dim > const &tria)
Definition adaptive_mesh_refinement.h:128
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:81
Definition adaptive_mesh_refinement.h:44