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
24#include <memory>
25
27
28template <unsigned int dim, unsigned int degree, typename number>
30{
31public:
38 , num_quad_points(SystemWide<dim, degree>::quadrature.size())
40 solve_context->get_user_inputs().spatial_discretization.max_refinement)
42 solve_context->get_user_inputs().spatial_discretization.min_refinement)
44 {
45 fe_values_flags.fill(dealii::UpdateFlags::update_default);
46 std::map<std::string, Types::Index> field_indices =
47 field_index_map(solve_context->get_field_attributes());
48 for (const auto &[name, field_criterion] :
49 solve_context->get_user_inputs().spatial_discretization.refinement_criteria)
50 {
51 // Grab the index and field type
52 const unsigned field_index = field_indices.at(name);
53 const TensorRank rank =
54 solve_context->get_field_attributes().at(field_index).field_type;
55
57 {
58 fe_values_flags[int(rank)] |= dealii::UpdateFlags::update_values;
59 }
60 else if (field_criterion.criterion & RefinementFlags::Gradient)
61 {
62 fe_values_flags[int(rank)] |= dealii::UpdateFlags::update_gradients;
63 }
64 }
65 }
66
70 ~RefinementManager() = default;
71
77 RefinementManager(const RefinementManager &grid_refiner) = delete;
78
85 operator=(const RefinementManager &grid_refiner) = delete;
86
92 RefinementManager(RefinementManager &&grid_refiner) noexcept = delete;
93
100 operator=(RefinementManager &&grid_refiner) noexcept = delete;
101
108 void
110 std::vector<std::shared_ptr<SolverBase<dim, degree, number>>> &solvers)
111 {
112 // Return early if adaptive meshing is disabled
113 if (!solve_context->get_user_inputs().spatial_discretization.has_adaptivity)
114 {
115 return;
116 }
117
118 // Mark cells and refine until no more cells are marked
120 bool first_iteration = true;
121 while (
122 dealii::Utilities::MPI::logical_or(mark_cells_for_refinement(), MPI_COMM_WORLD) ||
124 {
125 first_iteration = false;
126 refine_grid(solvers);
127 }
128
129 // Update anything affected by the grid change:
130 // Recalculate InvM since the grid has changed
131 solve_context->get_invm_manager().reinit(solve_context->get_dof_manager(),
132 solve_context->get_constraint_manager());
133 solve_context->get_invm_manager().compute_invm();
134 // Update the ghosts
135 Timer::start_section("Update ghosts");
136 for (auto &solver : solvers)
137 {
138 solver->update_ghosts();
139 }
140 Timer::end_section("Update ghosts");
141 }
142
143 void
145 {
146 marker_functions.push_back(marker);
147 }
148
149 void
151 {
152 marker_functions.clear();
153 }
154
155 const std::vector<std::shared_ptr<const CellMarkerBase<dim>>> &
157 {
158 return marker_functions;
159 }
160
167 void
169 std::vector<std::shared_ptr<SolverBase<dim, degree, number>>> &solvers)
170 {
172 solve_context->get_user_inputs().spatial_discretization;
173 const DoFManager<dim, degree> &dof_manager = solve_context->get_dof_manager();
174
175 dealii::types::global_dof_index old_dofs = dof_manager.get_total_dofs();
176 dealii::types::global_dof_index new_dofs = 0;
177 for (unsigned int remesh_index = 0; remesh_index < (space_parameters.max_refinement -
178 space_parameters.min_refinement);
179 remesh_index++)
180 {
181 // Perform grid refinement
182 ConditionalOStreams::pout_base() << "performing grid refinement...\n"
183 << std::flush;
184 do_adaptive_refinement(solvers);
185 // Update the ghosts
186 Timer::start_section("Update ghosts");
187 for (auto &solver : solvers)
188 {
189 solver->update_ghosts();
190 }
191 Timer::end_section("Update ghosts");
192
193 // Recalculate the total DoFs
194 new_dofs = dof_manager.get_total_dofs();
195
196 // Check for convergence
197 if (old_dofs == new_dofs)
198 {
199 break;
200 }
202 }
203 }
204
205private:
209 void
211 {
212 // Create the an object for the refinement criterion at each of the quad points. This
213 // will either contain the value for scalar fields, the magnitude for vector fields,
214 // or the magnitude of the gradient for both of the fields.
215 std::vector<number> values(num_quad_points, 0.0);
216
217 // Clear user flags
218 solve_context->get_triangulation_manager().clear_user_flags();
219
220 // Loop over the cells provided by the triangulation
221 for (const auto &cell : solve_context->get_triangulation_manager()
222 .get_triangulation()
223 .active_cell_iterators())
224 {
225 if (cell->is_locally_owned())
226 {
227 // Whether we should refine the cell
228 bool should_refine = false;
229
230 // TODO (landinjm): We can probably avoid checking some of the neighboring
231 // cells when coarsening them
232 std::map<std::string, Types::Index> field_indices =
233 field_index_map(solve_context->get_field_attributes());
234 for (const auto &[name, field_criterion] :
235 solve_context->get_user_inputs()
236 .spatial_discretization.refinement_criteria)
237 {
238 // Grab the index
239 const unsigned int index = field_indices.at(name);
240
241 // Grab the field type
243 solve_context->get_field_attributes().at(index).field_type;
244
245 // Grab the DoFHandler iterator
246 const auto dof_iterator = cell->as_dof_handler_iterator(
247 solve_context->get_dof_manager().get_field_dof_handler(index));
248
249 // Reinit the cell
250 dealii::FEValues<dim> fe_values(
254
255 fe_values.reinit(dof_iterator);
256 const auto &solution_vector =
257 solve_context->get_solution_indexer().get_solution_vector(index);
258
260 {
262 {
263 // Get the values for a scalar field
264 fe_values.get_function_values(solution_vector, values);
265 }
266 else
267 {
268 // Get the magnitude of the value for vector fields
269 // TODO (landinjm): Should be zeroing this out?
270 std::vector<dealii::Vector<number>> vector_values(
272 dealii::Vector<number>(dim));
273 fe_values.get_function_values(solution_vector, vector_values);
274 for (unsigned int q_point = 0; q_point < num_quad_points;
275 ++q_point)
276 {
277 values[q_point] = vector_values[q_point].l2_norm();
278 }
279 }
280
281 // Check if any of the quadrature points meet the refinement criterion
282 for (unsigned int q_point = 0; q_point < num_quad_points; ++q_point)
283 {
284 if (field_criterion.value_in_open_range(values[q_point]))
285 {
286 should_refine = true;
287 break;
288 }
289 }
290
291 // Exit if we've already determined that the cell has to be refined
292 if (should_refine)
293 {
294 break;
295 }
296 }
298 {
300 {
301 // Get the magnitude of the gradient for a scalar field
302 // TODO (landinjm): Should be zeroing this out?
303 std::vector<dealii::Tensor<1, dim, number>> scalar_gradients(
305 fe_values.get_function_gradients(solution_vector,
307 for (unsigned int q_point = 0; q_point < num_quad_points;
308 ++q_point)
309 {
310 values[q_point] = scalar_gradients[q_point].norm();
311 }
312 }
313 else
314 {
315 // TODO (landinjm): Should be zeroing this out?
316 std::vector<std::vector<dealii::Tensor<1, dim, number>>>
318 std::vector<dealii::Tensor<1, dim, number>>(
319 dim));
320 fe_values.get_function_gradients(solution_vector,
322 for (unsigned int q_point = 0; q_point < num_quad_points;
323 ++q_point)
324 {
325 dealii::Vector<number> vector_gradient_component_magnitude(
326 dim);
327 for (unsigned int dimension = 0; dimension < dim; dimension++)
328 {
330 vector_gradients[q_point][dimension].norm();
331 }
332 values[q_point] =
334 }
335 }
336
337 // Check if any of the quadrature points meet the refinement criterion
338 for (unsigned int q_point = 0; q_point < num_quad_points; ++q_point)
339 {
340 if (field_criterion.gradient_magnitude_above_threshold(
341 values[q_point]))
342 {
343 should_refine = true;
344 break;
345 }
346 }
347
348 // Exit if we've already determined that the cell has to be refined
349 if (should_refine)
350 {
351 break;
352 }
353 }
354 }
355
356 Assert(cell->level() > 0,
357 dealii::ExcMessage("Cell refinement level is less than one, which "
358 "will lead to underflow."));
359 const auto cell_refinement = static_cast<unsigned int>(cell->level());
361 {
362 cell->set_user_flag();
363 cell->clear_coarsen_flag();
364 cell->set_refine_flag();
365 }
366 if (should_refine)
367 {
368 cell->set_user_flag();
369 cell->clear_coarsen_flag();
370 }
372 !cell->user_flag_set())
373 {
374 cell->set_coarsen_flag();
375 }
376 }
377 }
378 }
379
386 bool
388 {
389 bool any_cell_marked = false;
390 for (const auto &cell : solve_context->get_triangulation_manager()
391 .get_triangulation()
392 .active_cell_iterators())
393 {
394 if (cell->is_locally_owned())
395 {
396 const unsigned int cell_refinement = cell->level();
397 if (std::any_of(
398 marker_functions.begin(),
399 marker_functions.end(),
400 [&](const std::shared_ptr<const CellMarkerBase<dim>> &marker_function)
401 {
402 return marker_function->flag(*cell,
403 solve_context->get_simulation_timer());
404 }))
405 {
406 cell->set_user_flag();
407 cell->clear_coarsen_flag();
409 {
410 cell->set_refine_flag();
411 any_cell_marked = true;
412 }
413 }
414 }
415 }
416 return any_cell_marked;
417 }
418
422 void
423 refine_grid(std::vector<std::shared_ptr<SolverBase<dim, degree, number>>> &solvers)
424 {
425 TriangulationManager<dim> &triangulation_manager =
426 solve_context->get_triangulation_manager();
427 DoFManager<dim, degree> &dof_manager = solve_context->get_dof_manager();
428 ConstraintManager<dim, degree, number> &constraint_manager =
429 solve_context->get_constraint_manager();
430
431 // Update ghosts of all fields.
432 for (auto solver : solvers)
433 {
434 solver->update_ghosts();
435 }
436
437 // Prepare for grid refinement
438 triangulation_manager.prepare_for_grid_refinement();
439 for (auto &solver : solvers)
440 {
441 solver->prepare_for_solution_transfer();
442 }
443
444 // Execute grid refinement
445 triangulation_manager.execute_grid_refinement();
446
447 // Redistribute DoFs and reinit the solvers
448 triangulation_manager.reinit();
449 dof_manager.reinit(triangulation_manager);
450 constraint_manager.reinit(solve_context->get_field_attributes());
451
452 // Reinit solutions, apply constraints, then solution transfer
453 for (auto &solver : solvers)
454 {
455 solver->reinit();
456 solver->execute_solution_transfer();
457 }
458 }
459
464
470 std::array<dealii::UpdateFlags, 2> fe_values_flags;
471
475 unsigned int num_quad_points = 0;
476
480 unsigned int max_refinement = 0;
481
485 unsigned int min_refinement = 0;
486
490 std::list<std::shared_ptr<const CellMarkerBase<dim>>> marker_functions;
491};
492
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:43
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:98
Class that manages the deal.II DoFHandlers.
Definition dof_manager.h:27
dealii::types::global_dof_index get_total_dofs() const
Get the total DoFs excluding multigrid DoFs.
Definition dof_manager.h:90
void reinit(const TriangulationManager< dim > &triangulation_handler)
Reinitialize the DoFHandlers.
Definition dof_manager.cc:59
Definition refinement_manager.h:30
unsigned int min_refinement
Minimum global refinement level.
Definition refinement_manager.h:485
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:387
unsigned int max_refinement
Maximum global refinement level.
Definition refinement_manager.h:480
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:470
void mark_cells_for_refinement_and_coarsening()
Mark cells for refinement and coarsening.
Definition refinement_manager.h:210
RefinementManager(SolveContext< dim, degree, number > &_solve_context)
Constructor. Init the flags for refinement.
Definition refinement_manager.h:35
RefinementManager(RefinementManager &&grid_refiner) noexcept=delete
Move constructor.
void clear_refinement_markers()
Definition refinement_manager.h:150
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:144
const std::vector< std::shared_ptr< const CellMarkerBase< dim > > > & get_refinement_markers() const
Definition refinement_manager.h:156
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:168
void do_adaptive_refinement(std::vector< std::shared_ptr< SolverBase< dim, degree, number > > > &solvers)
Do the adaptive refinement.
Definition refinement_manager.h:109
unsigned int num_quad_points
Number of quadrature points.
Definition refinement_manager.h:475
SolveContext< dim, degree, number > * solve_context
Grid refinement context.
Definition refinement_manager.h:463
void refine_grid(std::vector< std::shared_ptr< SolverBase< dim, degree, number > > > &solvers)
Refine the grid once.
Definition refinement_manager.h:423
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:490
This class provides context for a solver with ptrs to all the relevant dependencies.
Definition solve_context.h:36
Definition solver_base.h:32
This is the main class that handles the construction and solving of user-specified PDEs.
Definition system_wide.h:24
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:230
TensorRank
Tensor rank of the field.
Definition type_enums.h:30
@ Scalar
Definition type_enums.h:32