6#include <deal.II/base/bounding_box.h>
7#include <deal.II/fe/fe_values.h>
22#include <prismspf/config.h>
30template <
unsigned int dim,
unsigned int degree,
typename number>
48 std::map<std::string, Types::Index> field_indices =
50 for (
const auto &[name, field_criterion] :
51 solve_context->get_user_inputs().spatial_discretization.refinement_criteria)
54 const unsigned field_index = field_indices.at(name);
56 solve_context->get_field_attributes().at(field_index).field_type;
115 if (!
solve_context->get_user_inputs().spatial_discretization.has_adaptivity)
122 bool first_iteration =
true;
127 first_iteration =
false;
138 for (
auto &solver : solvers)
140 solver->update_ghosts();
157 const std::vector<std::shared_ptr<const CellMarkerBase<dim>>> &
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 -
189 for (
auto &solver : solvers)
191 solver->update_ghosts();
199 if (old_dofs == new_dofs)
220 solve_context->get_triangulation_manager().clear_user_flags();
223 for (
const auto &cell :
solve_context->get_triangulation_manager()
225 .active_cell_iterators())
227 if (cell->is_locally_owned())
230 bool should_refine =
false;
234 std::map<std::string, Types::Index> field_indices =
236 for (
const auto &[name, field_criterion] :
238 .spatial_discretization.refinement_criteria)
241 const unsigned int index = field_indices.at(name);
248 const auto dof_iterator = cell->as_dof_handler_iterator(
249 solve_context->get_dof_manager().get_field_dof_handler(index));
252 dealii::FEValues<dim> fe_values(
257 fe_values.reinit(dof_iterator);
258 const auto &solution_vector =
259 solve_context->get_solution_indexer().get_solution_vector(index);
266 fe_values.get_function_values(solution_vector, values);
272 std::vector<dealii::Vector<number>> vector_values(
274 dealii::Vector<number>(dim));
275 fe_values.get_function_values(solution_vector, vector_values);
279 values[q_point] = vector_values[q_point].l2_norm();
286 if (field_criterion.value_in_open_range(values[q_point]))
288 should_refine =
true;
305 std::vector<dealii::Tensor<1, dim, number>> scalar_gradients(
307 fe_values.get_function_gradients(solution_vector,
312 values[q_point] = scalar_gradients[q_point].norm();
318 std::vector<std::vector<dealii::Tensor<1, dim, number>>>
320 std::vector<dealii::Tensor<1, dim, number>>(
322 fe_values.get_function_gradients(solution_vector,
327 dealii::Vector<number> vector_gradient_component_magnitude(
329 for (
unsigned int dimension = 0; dimension < dim; dimension++)
331 vector_gradient_component_magnitude[dimension] =
332 vector_gradients[q_point][dimension].norm();
335 vector_gradient_component_magnitude.l2_norm();
342 if (field_criterion.gradient_magnitude_above_threshold(
345 should_refine =
true;
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());
364 cell->set_user_flag();
365 cell->clear_coarsen_flag();
366 cell->set_refine_flag();
370 cell->set_user_flag();
371 cell->clear_coarsen_flag();
374 !cell->user_flag_set())
376 cell->set_coarsen_flag();
391 bool any_cell_marked =
false;
392 for (
const auto &cell :
solve_context->get_triangulation_manager()
394 .active_cell_iterators())
396 if (cell->is_locally_owned())
398 const unsigned int cell_refinement = cell->level();
404 return marker_function->flag(*cell,
405 solve_context->get_simulation_timer());
408 cell->set_user_flag();
409 cell->clear_coarsen_flag();
412 cell->set_refine_flag();
413 any_cell_marked =
true;
418 return any_cell_marked;
436 for (
auto solver : solvers)
438 solver->update_ghosts();
443 for (
auto &solver : solvers)
445 solver->prepare_for_solution_transfer();
452 triangulation_manager.
reinit();
453 dof_manager.
reinit(triangulation_manager);
455 matrix_free_manager.
reinit(dof_manager, constraint_manager);
458 for (
auto &solver : solvers)
461 solver->execute_solution_transfer();
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