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/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
40namespace ExaDG
41{
42struct AdaptiveMeshRefinementData
43{
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)
51 {
52 }
53
54 void
55 print(dealii::ConditionalOStream const & pcout) const
56 {
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);
64 }
65
66 unsigned int trigger_every_n_time_steps;
67
68 int maximum_refinement_level;
69 int minimum_refinement_level;
70 bool preserve_boundary_cells;
71
72 double fraction_of_cells_to_be_refined;
73 double fraction_of_cells_to_be_coarsened;
74};
75
79inline bool
80trigger_coarsening_and_refinement_now(unsigned int const trigger_every_n_time_steps,
81 unsigned int const time_step_number)
82{
83 return ((time_step_number == 0) or not(time_step_number % trigger_every_n_time_steps));
84}
85
91template<int dim>
92void
93limit_coarsening_and_refinement(dealii::Triangulation<dim> & tria,
94 AdaptiveMeshRefinementData const & amr_data)
95{
96 for(auto & cell : tria.active_cell_iterators())
97 {
98 if(cell->is_locally_owned())
99 {
100 // Clear refinement flags on maximum refinement level.
101 if(cell->level() == amr_data.maximum_refinement_level)
102 {
103 cell->clear_refine_flag();
104 }
105
106 // Clear coarsening flags on minimum refinement level.
107 if(cell->level() <= amr_data.minimum_refinement_level)
108 {
109 cell->clear_coarsen_flag();
110 }
111
112 // Do not coarsen/refine cells along boundary if requested
113 if(amr_data.preserve_boundary_cells and cell->at_boundary())
114 {
115 cell->clear_refine_flag();
116 cell->clear_coarsen_flag();
117 }
118 }
119 }
120}
121
125template<int dim>
126bool
127any_cells_flagged_for_coarsening_or_refinement(dealii::Triangulation<dim> const & tria)
128{
129 bool any_flag_set = false;
130
131 for(auto & cell : tria.active_cell_iterators())
132 {
133 if(cell->is_locally_owned())
134 {
135 if(cell->refine_flag_set() or cell->coarsen_flag_set())
136 {
137 any_flag_set = true;
138 break;
139 }
140 }
141 }
142
143 any_flag_set = dealii::Utilities::MPI::logical_or(any_flag_set, tria.get_communicator());
144
145 return any_flag_set;
146}
147
148template<int dim, typename Number, typename VectorType>
149void
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)
157{
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();
165
166 dealii::QGauss<dim - 1> face_quadrature(n_face_quadrature_points);
167
168 dealii::Vector<float> estimated_error_per_cell(tria.n_active_cells());
169
170 dealii::KellyErrorEstimator<dim>::estimate(mapping,
171 dof_handler,
172 face_quadrature,
173 {}, // Neumann BC
174 locally_relevant_solution,
175 estimated_error_per_cell);
176
177 dealii::parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number(
178 tria,
179 estimated_error_per_cell,
180 amr_data.fraction_of_cells_to_be_refined,
181 amr_data.fraction_of_cells_to_be_coarsened);
182}
183
184} // namespace ExaDG
185
186#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: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