4#ifndef triangulation_handler_h
5#define triangulation_handler_h
7#include <deal.II/distributed/tria.h>
8#include <deal.II/grid/tria.h>
10#include <prismspf/config.h>
11#include <prismspf/core/solution_handler.h>
12#include <prismspf/user_inputs/user_input_parameters.h>
17PRISMS_PF_BEGIN_NAMESPACE
27 std::conditional_t<dim == 1,
28 dealii::Triangulation<dim>,
29 dealii::parallel::distributed::Triangulation<dim>>;
39 [[nodiscard]]
const Triangulation &
45 [[nodiscard]]
const std::vector<std::shared_ptr<const dealii::Triangulation<dim>>> &
52 [[nodiscard]]
const dealii::Triangulation<dim> &
58 [[nodiscard]]
unsigned int
64 [[nodiscard]]
unsigned int
70 [[nodiscard]]
unsigned int
101 mark_boundaries()
const;
112 const userInputParameters<dim> *user_inputs;
117 std::shared_ptr<Triangulation> triangulation;
126 std::vector<std::shared_ptr<const dealii::Triangulation<dim>>> coarsened_triangulations;
131 bool has_multigrid =
false;
136 unsigned int min_level = 0;
141 unsigned int max_level = 0;
144PRISMS_PF_END_NAMESPACE
This class handlers the generation and manipulation of triangulations.
Definition triangulation_handler.h:24
const std::vector< std::shared_ptr< const dealii::Triangulation< dim > > > & get_mg_triangulation() const
Getter function for the multigrid triangulation (constant reference).
Definition triangulation_handler.cc:54
void export_triangulation_as_vtk(const std::string &filename) const
Export triangulation to vtk. This is done for debugging purposes when dealing with unusual meshes (e....
Definition triangulation_handler.cc:187
unsigned int get_mg_min_level() const
Return the maximum multigrid level.
Definition triangulation_handler.cc:83
void generate_mesh()
Generate mesh based on the inputs provided by the user.
Definition triangulation_handler.cc:101
unsigned int get_n_global_levels() const
Return the global maximum level of the triangulation.
Definition triangulation_handler.cc:75
const Triangulation & get_triangulation() const
Getter function for triangulation (constant reference).
Definition triangulation_handler.cc:46
unsigned int get_mg_max_level() const
Return the minimum multigrid level.
Definition triangulation_handler.cc:92
void adaptively_refine_mesh(solutionHandler< dim > &solution_handler)
Adaptively refine the mesh based on the inputs provided by the user.
Definition triangulation_handler.cc:178