![]() |
ryujin 2.1.1 revision dbf0e3ba7acdb60b6d558e4257815df4a8f8daf9
|
#include <source/mesh_adaptor.h>
Public Member Functions | |
Constructor and setup | |
MeshAdaptor (const MPIEnsemble &mpi_ensemble, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem ¶bolic_system, const InitialPrecomputedVector &initial_precomputed, const ScalarVector &alpha, const std::string &subsection="/MeshAdaptor") | |
void | prepare (const Number t) |
void | analyze (const StateVector &state_vector, const Number t, unsigned int cycle) |
const auto & | need_mesh_adaptation () const |
void | mark_cells_for_coarsening_and_refinement (dealii::Triangulation< dim > &triangulation) const |
const auto & | indicators () const |
void | compute_smoothness_indicators (const StateVector &state_vector) const |
const auto & | smoothness_indicators () const |
Typedefs and constexpr constants | |
using | HyperbolicSystem = typename Description::HyperbolicSystem |
using | ParabolicSystem = typename Description::ParabolicSystem |
using | View = typename Description::template HyperbolicSystemView< dim, Number > |
using | StateVector = typename View::StateVector |
using | InitialPrecomputedVector = typename View::InitialPrecomputedVector |
using | ScalarVector = Vectors::ScalarVector< Number > |
static constexpr auto | problem_dimension = View::problem_dimension |
The MeshAdaptor class is responsible for performing global or local mesh adaptation.
Definition at line 106 of file mesh_adaptor.h.
using ryujin::MeshAdaptor< Description, dim, Number >::HyperbolicSystem = typename Description::HyperbolicSystem |
Definition at line 114 of file mesh_adaptor.h.
using ryujin::MeshAdaptor< Description, dim, Number >::ParabolicSystem = typename Description::ParabolicSystem |
Definition at line 115 of file mesh_adaptor.h.
using ryujin::MeshAdaptor< Description, dim, Number >::View = typename Description::template HyperbolicSystemView<dim, Number> |
Definition at line 117 of file mesh_adaptor.h.
using ryujin::MeshAdaptor< Description, dim, Number >::StateVector = typename View::StateVector |
Definition at line 122 of file mesh_adaptor.h.
using ryujin::MeshAdaptor< Description, dim, Number >::InitialPrecomputedVector = typename View::InitialPrecomputedVector |
Definition at line 123 of file mesh_adaptor.h.
using ryujin::MeshAdaptor< Description, dim, Number >::ScalarVector = Vectors::ScalarVector<Number> |
Definition at line 124 of file mesh_adaptor.h.
ryujin::MeshAdaptor< Description, dim, Number >::MeshAdaptor | ( | const MPIEnsemble & | mpi_ensemble, |
const OfflineData< dim, Number > & | offline_data, | ||
const HyperbolicSystem & | hyperbolic_system, | ||
const ParabolicSystem & | parabolic_system, | ||
const InitialPrecomputedVector & | initial_precomputed, | ||
const ScalarVector & | alpha, | ||
const std::string & | subsection = "/MeshAdaptor< Description, dim, Number >" |
||
) |
Constructor.
Definition at line 22 of file mesh_adaptor.template.h.
References ryujin::fixed_threshold, ryujin::simulation_cycle, and ryujin::smoothness_indicators.
void ryujin::MeshAdaptor< Description, dim, Number >::prepare | ( | const Number | t | ) |
Prepare temporary storage and clean up internal data for the analyze() facility.
Definition at line 166 of file mesh_adaptor.template.h.
void ryujin::MeshAdaptor< Description, dim, Number >::analyze | ( | const StateVector & | state_vector, |
const Number | t, | ||
unsigned int | cycle | ||
) |
Analyze the given StateVector with the configured adaptation strategy and time point selection strategy and decide whether a mesh adaptation cycle should be performed.
Definition at line 412 of file mesh_adaptor.template.h.
References ryujin::fixed_time_points, ryujin::global_refinement, ryujin::random_adaptation, ryujin::simulation_cycle, and ryujin::smoothness_indicators.
|
inline |
A boolean indicating whether we should perform a mesh adapation step in the current cycle. The analyze() method will set this boolean to true whenever the selected adaptation strategy advices to perform an adaptation cycle.
Definition at line 164 of file mesh_adaptor.h.
void ryujin::MeshAdaptor< Description, dim, Number >::mark_cells_for_coarsening_and_refinement | ( | dealii::Triangulation< dim > & | triangulation | ) | const |
Mark cells for coarsening and refinement with the configured mesh adaptation and marking strategies.
Definition at line 528 of file mesh_adaptor.template.h.
References ryujin::fixed_threshold, ryujin::global_refinement, ryujin::maximum, ryujin::minimum, ryujin::random_adaptation, and ryujin::smoothness_indicators.
|
inline |
The computed cell indicators.
Definition at line 176 of file mesh_adaptor.h.
void ryujin::MeshAdaptor< Description, dim, Number >::compute_smoothness_indicators | ( | const StateVector & | state_vector | ) | const |
Compute smoothness indicators. This function reinitializes and populates the smoothness_indicators() vector.
Definition at line 191 of file mesh_adaptor.template.h.
References ryujin::SelectedComponentsExtractor< Description, dim, Number >::extract(), RYUJIN_OMP_FOR, RYUJIN_PARALLEL_REGION_BEGIN, and RYUJIN_PARALLEL_REGION_END.
|
inline |
The computed smoothness indicators. The vector is only valid if the "smoothness indicator" refinement strategy has been selected.
Definition at line 188 of file mesh_adaptor.h.
|
staticconstexpr |
Definition at line 120 of file mesh_adaptor.h.