PRISMS-PF Manual
Loading...
Searching...
No Matches
refinement_manager.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/bounding_box.h>
7#include <deal.II/fe/fe_values.h>
8
14#include <prismspf/core/timer.h>
16
19
21
22#include <prismspf/config.h>
23
25
26#include <memory>
27
29
30template <unsigned int dim, unsigned int degree, typename number>
32{
33public:
38 : solve_context(&_solve_context)
40 , num_quad_points(SystemWide<dim, degree>::quadrature.size())
42 solve_context->get_user_inputs().spatial_discretization.max_refinement)
44 solve_context->get_user_inputs().spatial_discretization.min_refinement)
46 {
47 fe_values_flags.fill(dealii::UpdateFlags::update_default);
48 std::map<std::string, Types::Index> field_indices =
49 field_index_map(solve_context->get_field_attributes());
50 for (const auto &[name, field_criterion] :
51 solve_context->get_user_inputs().spatial_discretization.refinement_criteria)
52 {
53 // Grab the index and field type
54 const unsigned field_index = field_indices.at(name);
55 const TensorRank rank =
56 solve_context->get_field_attributes().at(field_index).field_type;
57
58 if (field_criterion.criterion & RefinementFlags::Value)
59 {
60 fe_values_flags[int(rank)] |= dealii::UpdateFlags::update_values;
61 }
62 else if (field_criterion.criterion & RefinementFlags::Gradient)
63 {
64 fe_values_flags[int(rank)] |= dealii::UpdateFlags::update_gradients;
65 }
66 }
67 }
68
72 ~RefinementManager() = default;
73
79 RefinementManager(const RefinementManager &grid_refiner) = delete;
80
87 operator=(const RefinementManager &grid_refiner) = delete;
88
94 RefinementManager(RefinementManager &&grid_refiner) noexcept = delete;
95
102 operator=(RefinementManager &&grid_refiner) noexcept = delete;
103
110 void
112 std::vector<std::shared_ptr<SolverBase<dim, degree, number>>> &solvers)
113 {
114 // Return early if adaptive meshing is disabled
115 if (!solve_context->get_user_inputs().spatial_discretization.has_adaptivity)
116 {
117 return;
118 }
119
120 // Mark cells and refine until no more cells are marked
122 bool first_iteration = true;
123 while (
124 dealii::Utilities::MPI::logical_or(mark_cells_for_refinement(), MPI_COMM_WORLD) ||
125 first_iteration)
126 {
127 first_iteration = false;
128 refine_grid(solvers);
129 }
130
131 // Update anything affected by the grid change:
132 // Recalculate InvM since the grid has changed
133 solve_context->get_invm_manager().reinit(solve_context->get_dof_manager(),
134 solve_context->get_constraint_manager());
135 solve_context->get_invm_manager().compute_invm();
136 // Update the ghosts
137 Timer::start_section("Update ghosts");
138 for (auto &solver : solvers)
139 {
140 solver->update_ghosts();
141 }
142 Timer::end_section("Update ghosts");
143 }
144
145 void
146 add_refinement_marker(std::shared_ptr<const CellMarkerBase<dim>> marker)
147 {
148 marker_functions.push_back(marker);
149 }
150
151 void
153 {
154 marker_functions.clear();
155 }
156
157 const std::vector<std::shared_ptr<const CellMarkerBase<dim>>> &
159 {
160 return marker_functions;
161 }
162
169 void
171 std::vector<std::shared_ptr<SolverBase<dim, degree, number>>> &solvers)
172 {
173 const SpatialDiscretization<dim> &space_parameters =
174 solve_context->get_user_inputs().spatial_discretization;
175 const DoFManager<dim, degree> &dof_manager = solve_context->get_dof_manager();
176
177 dealii::types::global_dof_index old_dofs = dof_manager.get_total_dofs();
178 dealii::types::global_dof_index new_dofs = 0;
179 for (unsigned int remesh_index = 0; remesh_index < (space_parameters.max_refinement -
180 space_parameters.min_refinement);
181 remesh_index++)
182 {
183 // Perform grid refinement
184 ConditionalOStreams::pout_base() << "performing grid refinement...\n"
185 << std::flush;
186 do_adaptive_refinement(solvers);
187 // Update the ghosts
188 Timer::start_section("Update ghosts");
189 for (auto &solver : solvers)
190 {
191 solver->update_ghosts();
192 }
193 Timer::end_section("Update ghosts");
194
195 // Recalculate the total DoFs
196 new_dofs = dof_manager.get_total_dofs();
197
198 // Check for convergence
199 if (old_dofs == new_dofs)
200 {
201 break;
202 }
203 old_dofs = new_dofs;
204 }
205 }
206
207private:
211 void
213 {
214 // Create the an object for the refinement criterion at each of the quad points. This
215 // will either contain the value for scalar fields, the magnitude for vector fields,
216 // or the magnitude of the gradient for both of the fields.
217 std::vector<number> values(num_quad_points, 0.0);
218
219 // Clear user flags
220 solve_context->get_triangulation_manager().clear_user_flags();
221
222 // Loop over the cells provided by the triangulation
223 for (const auto &cell : solve_context->get_triangulation_manager()
224 .get_triangulation()
225 .active_cell_iterators())
226 {
227 if (cell->is_locally_owned())
228 {
229 // Whether we should refine the cell
230 bool should_refine = false;
231
232 // TODO (landinjm): We can probably avoid checking some of the neighboring
233 // cells when coarsening them
234 std::map<std::string, Types::Index> field_indices =
235 field_index_map(solve_context->get_field_attributes());
236 for (const auto &[name, field_criterion] :
237 solve_context->get_user_inputs()
238 .spatial_discretization.refinement_criteria)
239 {
240 // Grab the index
241 const unsigned int index = field_indices.at(name);
242
243 // Grab the field type
244 const TensorRank local_field_type =
245 solve_context->get_field_attributes().at(index).field_type;
246
247 // Grab the DoFHandler iterator
248 const auto dof_iterator = cell->as_dof_handler_iterator(
249 solve_context->get_dof_manager().get_field_dof_handler(index));
250
251 // Reinit the cell
252 dealii::FEValues<dim> fe_values(
253 SystemWide<dim, degree>::fe_systems[int(local_field_type)],
255 fe_values_flags.at(int(local_field_type)));
256
257 fe_values.reinit(dof_iterator);
258 const auto &solution_vector =
259 solve_context->get_solution_indexer().get_solution_vector(index);
260
261 if (field_criterion.criterion & RefinementFlags::Value)
262 {
263 if (local_field_type == TensorRank::Scalar)
264 {
265 // Get the values for a scalar field
266 fe_values.get_function_values(solution_vector, values);
267 }
268 else
269 {
270 // Get the magnitude of the value for vector fields
271 // TODO (landinjm): Should be zeroing this out?
272 std::vector<dealii::Vector<number>> vector_values(
274 dealii::Vector<number>(dim));
275 fe_values.get_function_values(solution_vector, vector_values);
276 for (unsigned int q_point = 0; q_point < num_quad_points;
277 ++q_point)
278 {
279 values[q_point] = vector_values[q_point].l2_norm();
280 }
281 }
282
283 // Check if any of the quadrature points meet the refinement criterion
284 for (unsigned int q_point = 0; q_point < num_quad_points; ++q_point)
285 {
286 if (field_criterion.value_in_open_range(values[q_point]))
287 {
288 should_refine = true;
289 break;
290 }
291 }
292
293 // Exit if we've already determined that the cell has to be refined
294 if (should_refine)
295 {
296 break;
297 }
298 }
299 if (field_criterion.criterion & RefinementFlags::Gradient)
300 {
301 if (local_field_type == TensorRank::Scalar)
302 {
303 // Get the magnitude of the gradient for a scalar field
304 // TODO (landinjm): Should be zeroing this out?
305 std::vector<dealii::Tensor<1, dim, number>> scalar_gradients(
307 fe_values.get_function_gradients(solution_vector,
308 scalar_gradients);
309 for (unsigned int q_point = 0; q_point < num_quad_points;
310 ++q_point)
311 {
312 values[q_point] = scalar_gradients[q_point].norm();
313 }
314 }
315 else
316 {
317 // TODO (landinjm): Should be zeroing this out?
318 std::vector<std::vector<dealii::Tensor<1, dim, number>>>
319 vector_gradients(num_quad_points,
320 std::vector<dealii::Tensor<1, dim, number>>(
321 dim));
322 fe_values.get_function_gradients(solution_vector,
323 vector_gradients);
324 for (unsigned int q_point = 0; q_point < num_quad_points;
325 ++q_point)
326 {
327 dealii::Vector<number> vector_gradient_component_magnitude(
328 dim);
329 for (unsigned int dimension = 0; dimension < dim; dimension++)
330 {
331 vector_gradient_component_magnitude[dimension] =
332 vector_gradients[q_point][dimension].norm();
333 }
334 values[q_point] =
335 vector_gradient_component_magnitude.l2_norm();
336 }
337 }
338
339 // Check if any of the quadrature points meet the refinement criterion
340 for (unsigned int q_point = 0; q_point < num_quad_points; ++q_point)
341 {
342 if (field_criterion.gradient_magnitude_above_threshold(
343 values[q_point]))
344 {
345 should_refine = true;
346 break;
347 }
348 }
349
350 // Exit if we've already determined that the cell has to be refined
351 if (should_refine)
352 {
353 break;
354 }
355 }
356 }
357
358 Assert(cell->level() > 0,
359 dealii::ExcMessage("Cell refinement level is less than one, which "
360 "will lead to underflow."));
361 const auto cell_refinement = static_cast<unsigned int>(cell->level());
362 if (should_refine && cell_refinement < max_refinement)
363 {
364 cell->set_user_flag();
365 cell->clear_coarsen_flag();
366 cell->set_refine_flag();
367 }
368 if (should_refine)
369 {
370 cell->set_user_flag();
371 cell->clear_coarsen_flag();
372 }
373 if (!should_refine && cell_refinement > min_refinement &&
374 !cell->user_flag_set())
375 {
376 cell->set_coarsen_flag();
377 }
378 }
379 }
380 }
381
388 bool
390 {
391 bool any_cell_marked = false;
392 for (const auto &cell : solve_context->get_triangulation_manager()
393 .get_triangulation()
394 .active_cell_iterators())
395 {
396 if (cell->is_locally_owned())
397 {
398 const unsigned int cell_refinement = cell->level();
399 if (std::any_of(
400 marker_functions.begin(),
401 marker_functions.end(),
402 [&](const std::shared_ptr<const CellMarkerBase<dim>> &marker_function)
403 {
404 return marker_function->flag(*cell,
405 solve_context->get_simulation_timer());
406 }))
407 {
408 cell->set_user_flag();
409 cell->clear_coarsen_flag();
410 if (cell_refinement < max_refinement)
411 {
412 cell->set_refine_flag();
413 any_cell_marked = true;
414 }
415 }
416 }
417 }
418 return any_cell_marked;
419 }
420
424 void
425 refine_grid(std::vector<std::shared_ptr<SolverBase<dim, degree, number>>> &solvers)
426 {
427 TriangulationManager<dim> &triangulation_manager =
428 solve_context->get_triangulation_manager();
429 DoFManager<dim, degree> &dof_manager = solve_context->get_dof_manager();
430 ConstraintManager<dim, degree, number> &constraint_manager =
431 solve_context->get_constraint_manager();
432 MatrixFreeManager<dim, number> &matrix_free_manager =
433 solve_context->get_matrix_free_manager();
434
435 // Update ghosts of all fields.
436 for (auto solver : solvers)
437 {
438 solver->update_ghosts();
439 }
440
441 // Prepare for grid refinement
442 triangulation_manager.prepare_for_grid_refinement();
443 for (auto &solver : solvers)
444 {
445 solver->prepare_for_solution_transfer();
446 }
447
448 // Execute grid refinement
449 triangulation_manager.execute_grid_refinement();
450
451 // Redistribute DoFs and reinit the solvers
452 triangulation_manager.reinit();
453 dof_manager.reinit(triangulation_manager);
454 constraint_manager.reinit(solve_context->get_field_attributes());
455 matrix_free_manager.reinit(dof_manager, constraint_manager);
456
457 // Reinit solutions, apply constraints, then solution transfer
458 for (auto &solver : solvers)
459 {
460 solver->reinit();
461 solver->execute_solution_transfer();
462 }
463 }
464
469
475 std::array<dealii::UpdateFlags, 2> fe_values_flags;
476
480 unsigned int num_quad_points = 0;
481
485 unsigned int max_refinement = 0;
486
490 unsigned int min_refinement = 0;
491
495 std::list<std::shared_ptr<const CellMarkerBase<dim>>> marker_functions;
496};
497
498PRISMS_PF_END_NAMESPACE
Base class for cell markers.
Definition cell_marker_base.h:25
static dealii::ConditionalOStream & pout_base()
Generic parallel output stream. Used for essential information in release and debug mode.
Definition conditional_ostreams.cc:44
The class handles the generation and application of boundary conditions based on the user-inputs.
Definition constraint_manager.h:50
void reinit(const std::vector< FieldAttributes > &field_attributes)
Make constraints based on the inputs of the constructor.
Definition constraint_manager.cc:112
Class that manages the deal.II DoFHandlers.
Definition dof_manager.h:25
dealii::types::global_dof_index get_total_dofs() const
Get the total DoFs excluding multigrid DoFs.
Definition dof_manager.cc:115
void reinit(const TriangulationManager< dim > &triangulation_manager)
Reinitialize the DoFHandlers.
Definition dof_manager.cc:33
Containers for matrix free objects.
Definition matrix_free_manager.h:28
void reinit(const DoFManager< dim, degree > &dof_manager, const ConstraintManager< dim, degree, number > &constraint_manager)
Reinit.
Definition matrix_free_manager.h:74
unsigned int min_refinement
Minimum global refinement level.
Definition refinement_manager.h:490
RefinementManager(const RefinementManager &grid_refiner)=delete
Copy constructor.
bool mark_cells_for_refinement()
Mark cells based on function. Note: cells are only marked for refinement but not coarsening.
Definition refinement_manager.h:389
unsigned int max_refinement
Maximum global refinement level.
Definition refinement_manager.h:485
std::array< dealii::UpdateFlags, 2 > fe_values_flags
Update flags for the FEValues determined by the grid refinement criterion. For now,...
Definition refinement_manager.h:475
void mark_cells_for_refinement_and_coarsening()
Mark cells for refinement and coarsening.
Definition refinement_manager.h:212
RefinementManager(SolveContext< dim, degree, number > &_solve_context)
Constructor. Init the flags for refinement.
Definition refinement_manager.h:37
RefinementManager(RefinementManager &&grid_refiner) noexcept=delete
Move constructor.
void clear_refinement_markers()
Definition refinement_manager.h:152
RefinementManager & operator=(RefinementManager &&grid_refiner) noexcept=delete
Move assignment.
~RefinementManager()=default
Destructor.
void add_refinement_marker(std::shared_ptr< const CellMarkerBase< dim > > marker)
Definition refinement_manager.h:146
const std::vector< std::shared_ptr< const CellMarkerBase< dim > > > & get_refinement_markers() const
Definition refinement_manager.h:158
void do_initial_refinement(std::vector< std::shared_ptr< SolverBase< dim, degree, number > > > &solvers)
Similar to do_adaptive_refinement but loops coarsening.
Definition refinement_manager.h:170
void do_adaptive_refinement(std::vector< std::shared_ptr< SolverBase< dim, degree, number > > > &solvers)
Do the adaptive refinement.
Definition refinement_manager.h:111
unsigned int num_quad_points
Number of quadrature points.
Definition refinement_manager.h:480
SolveContext< dim, degree, number > * solve_context
Grid refinement context.
Definition refinement_manager.h:468
void refine_grid(std::vector< std::shared_ptr< SolverBase< dim, degree, number > > > &solvers)
Refine the grid once.
Definition refinement_manager.h:425
RefinementManager & operator=(const RefinementManager &grid_refiner)=delete
Copy assignment.
std::list< std::shared_ptr< const CellMarkerBase< dim > > > marker_functions
Marker functions.
Definition refinement_manager.h:495
This class provides context for a solver with ptrs to all the relevant dependencies.
Definition solve_context.h:34
Definition solver_base.h:30
This is the main class that handles the construction and solving of user-specified PDEs.
Definition system_wide.h:24
static const std::array< const dealii::FESystem< dim >, 2 > fe_systems
Scalar and Vector FE systems.
Definition system_wide.h:29
static const dealii::QGaussLobatto< dim > quadrature
Quadrature rule.
Definition system_wide.h:41
static void start_section(const char *name)
Start a new timer section.
Definition timer.cc:116
static void end_section(const char *name)
End the timer section.
Definition timer.cc:127
This class handlers the generation and manipulation of triangulations.
Definition triangulation_manager.h:26
void reinit()
Reinitialize the triangulation handler. This is used for AMR with multigrid so the coarsened meshes c...
Definition triangulation_manager.cc:44
void execute_grid_refinement()
Execute grid refinement on the triangulation.
Definition triangulation_manager.h:109
void prepare_for_grid_refinement()
Prepare the triangulation for grid refinement.
Definition triangulation_manager.h:100
std::map< std::string, Types::Index > field_index_map(const std::vector< FieldAttributes > &fields)
Make a map that maps field names to field indices.
Definition field_attributes.h:64
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
@ Gradient
Use gradient of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:36
Definition conditional_ostreams.cc:20
Struct that holds spatial discretization parameters.
Definition spatial_discretization.h:231
unsigned int min_refinement
Definition spatial_discretization.h:327
unsigned int max_refinement
Definition spatial_discretization.h:324
TensorRank
Tensor rank of the field.
Definition type_enums.h:52
@ Scalar
Definition type_enums.h:54