Fluid structure interaction suite
ParsedTools::GridRefinement Class Reference

A wrapper for refinement strategies. More...

#include <grid_refinement.h>

Inheritance diagram for ParsedTools::GridRefinement:
[legend]

Public Member Functions

 GridRefinement (const std::string &section_name="", const unsigned int &n_refinement_cycles=1, const RefinementStrategy &strategy=RefinementStrategy::global, const std::string &estimator_type="kelly", const double &top_parameter=.3, const double &bottom_parameter=0.1, const unsigned int &max_cells=0, const int &min_level=0, const int &max_level=0, const std::map< std::string, std::function< void(Vector< float > &)>> &optional_estimators={}, const ComponentMask &component_mask=ComponentMask())
 Constructor. More...
 
template<int dim, int spacedim>
void mark_cells (const Vector< float > &criteria, Triangulation< dim, spacedim > &tria) const
 Mark cells a the triangulation for refinement or coarsening, according to the given strategy applied to the supplied vector representing local error criteria. More...
 
template<int dim, int spacedim>
void mark_cells (const Vector< float > &criteria, parallel::distributed::Triangulation< dim, spacedim > &tria) const
 Mark cells of a distribtued triangulation for refinement or coarsening, according to the given strategy applied to the supplied vector representing local error criteria. More...
 
template<int dim, int spacedim, typename VectorType >
void estimate_error (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Vector< float > &estimated_error_per_cell, const std::map< types::boundary_id, const Function< spacedim, typename VectorType::value_type > * > &neumann_bc={}, const Function< spacedim > *coefficients=nullptr, const typename KellyErrorEstimator< dim, spacedim >::Strategy strategy=KellyErrorEstimator< dim, spacedim >::cell_diameter_over_24) const
 Call the error estimator specified in the parameter file. More...
 
template<int dim, int spacedim, typename VectorType >
void estimate_error (const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Vector< float > &estimated_error_per_cell, const std::map< types::boundary_id, const Function< spacedim, typename VectorType::value_type > * > &neumann_bc={}, const Function< spacedim > *coefficients=nullptr, const typename KellyErrorEstimator< dim, spacedim >::Strategy strategy=KellyErrorEstimator< dim, spacedim >::cell_diameter_over_24) const
 Call the error estimator specified in the parameter file with a default mapping. More...
 
const RefinementStrategyget_strategy () const
 Get the current strategy object. More...
 
const unsigned int & get_n_refinement_cycles () const
 Get the total number of refinemt cycles. More...
 
std_cxx20::ranges::iota_view< unsigned int, unsigned int > get_refinement_cycles () const
 Get Return an object that can be thought of as an array containing all indices from zero to n_refinement_cycles. More...
 
template<int dim, int spacedim, typename VectorType , typename Tria >
void estimate_mark_refine (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Tria &tria) const
 Perform all the steps of the ESTIMATE-MARK-REFINE cycle. More...
 
template<int dim, int spacedim, typename VectorType , typename Tria >
void estimate_mark_refine (const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Tria &tria) const
 Perform all the steps of the ESTIMATE-MARK-REFINE cycle. More...
 
- Public Member Functions inherited from ParameterAcceptor
 ParameterAcceptor (const std::string &section_name="")
 
unsigned int get_acceptor_id () const
 
virtual ~ParameterAcceptor () override
 
virtual void declare_parameters (ParameterHandler &prm)
 
virtual void parse_parameters (ParameterHandler &prm)
 
std::string get_section_name () const
 
std::vector< std::string > get_section_path () const
 
void add_parameter (const std::string &entry, ParameterType &parameter, const std::string &documentation="", ParameterHandler &prm_=prm, const Patterns::PatternBase &pattern=*Patterns::Tools::Convert< ParameterType >::to_pattern())
 
void enter_subsection (const std::string &subsection)
 
void leave_subsection ()
 
void enter_my_subsection (ParameterHandler &prm)
 
void leave_my_subsection (ParameterHandler &prm)
 
void serialize (Archive &ar, const unsigned int version)
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void serialize (Archive &ar, const unsigned int version)
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Private Member Functions

template<int dim, int spacedim>
void limit_levels (Triangulation< dim, spacedim > &tria) const
 Make sure that the refinement level is kept between min_level and max_level. More...
 

Private Attributes

unsigned int n_refinement_cycles
 
RefinementStrategy strategy
 
double top_parameter
 
double bottom_parameter
 
unsigned int max_cells
 
int min_level
 
int max_level
 
std::string estimator_type
 
std::map< std::string, std::function< void(Vector< float > &)> > optional_estimators
 
ComponentMask component_mask
 
types::subdomain_id subdomain_id
 
types::material_id material_id
 

Additional Inherited Members

- Static Public Member Functions inherited from ParameterAcceptor
static void initialize (const std::string &filename="", const std::string &output_filename="", const ParameterHandler::OutputStyle output_style_for_output_filename=ParameterHandler::Short, ParameterHandler &prm=ParameterAcceptor::prm, const ParameterHandler::OutputStyle output_style_for_filename=ParameterHandler::DefaultStyle)
 
static void initialize (std::istream &input_stream, ParameterHandler &prm=ParameterAcceptor::prm)
 
static void clear ()
 
static void parse_all_parameters (ParameterHandler &prm=ParameterAcceptor::prm)
 
static void declare_all_parameters (ParameterHandler &prm=ParameterAcceptor::prm)
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
- Public Attributes inherited from ParameterAcceptor
boost::signals2::signal< void()> declare_parameters_call_back
 
boost::signals2::signal< void()> parse_parameters_call_back
 
- Static Public Attributes inherited from ParameterAcceptor
static ParameterHandler prm
 
- Protected Attributes inherited from ParameterAcceptor
const std::string section_name
 
std::vector< std::string > subsections
 

Detailed Description

A wrapper for refinement strategies.

This class implements a parametrized version of the ESTIMATE-MARK-REFINE steps in AFEM methods. In particular, it allows you to select between

and to use KellyErrorEstimator, or a custom estimator, to compute the error indicators in the first two cases.

In addition we have also some control parameters to guarantee that the number of cells is maintained under a maximum, and that the refinement level is kept between min_level and max_level.

Definition at line 68 of file grid_refinement.h.

Constructor & Destructor Documentation

◆ GridRefinement()

ParsedTools::GridRefinement::GridRefinement ( const std::string &  section_name = "",
const unsigned int &  n_refinement_cycles = 1,
const RefinementStrategy strategy = RefinementStrategy::global,
const std::string &  estimator_type = "kelly",
const double &  top_parameter = .3,
const double &  bottom_parameter = 0.1,
const unsigned int &  max_cells = 0,
const int &  min_level = 0,
const int &  max_level = 0,
const std::map< std::string, std::function< void(Vector< float > &)>  ,
optional_estimators = {},
const ComponentMask component_mask = ComponentMask() 
)

Member Function Documentation

◆ mark_cells() [1/2]

template<int dim, int spacedim>
void ParsedTools::GridRefinement::mark_cells ( const Vector< float > &  criteria,
Triangulation< dim, spacedim > &  tria 
) const

Mark cells a the triangulation for refinement or coarsening, according to the given strategy applied to the supplied vector representing local error criteria.

Cells are only marked for refinement or coarsening. No refinement is actually performed. You need to call Triangulation::execute_coarsening_and_refinement() yourself.

◆ mark_cells() [2/2]

template<int dim, int spacedim>
void ParsedTools::GridRefinement::mark_cells ( const Vector< float > &  criteria,
parallel::distributed::Triangulation< dim, spacedim > &  tria 
) const

Mark cells of a distribtued triangulation for refinement or coarsening, according to the given strategy applied to the supplied vector representing local error criteria.

If the criterion which is specified in the parameter file is not available, an exception is thrown.

Cells are only marked for refinement or coarsening. No refinement is actually performed. You need to call Triangulation::execute_coarsening_and_refinement() yourself.

◆ estimate_error() [1/2]

template<int dim, int spacedim, typename VectorType >
void ParsedTools::GridRefinement::estimate_error ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof_handler,
const VectorType &  solution,
Vector< float > &  estimated_error_per_cell,
const std::map< types::boundary_id, const Function< spacedim, typename VectorType::value_type > * > &  neumann_bc = {},
const Function< spacedim > *  coefficients = nullptr,
const typename KellyErrorEstimator< dim, spacedim >::Strategy  strategy = KellyErrorEstimator< dim, spacedim >::cell_diameter_over_24 
) const

Call the error estimator specified in the parameter file.

◆ estimate_error() [2/2]

template<int dim, int spacedim, typename VectorType >
void ParsedTools::GridRefinement::estimate_error ( const DoFHandler< dim, spacedim > &  dof_handler,
const VectorType &  solution,
Vector< float > &  estimated_error_per_cell,
const std::map< types::boundary_id, const Function< spacedim, typename VectorType::value_type > * > &  neumann_bc = {},
const Function< spacedim > *  coefficients = nullptr,
const typename KellyErrorEstimator< dim, spacedim >::Strategy  strategy = KellyErrorEstimator< dim, spacedim >::cell_diameter_over_24 
) const

Call the error estimator specified in the parameter file with a default mapping.

◆ get_strategy()

const RefinementStrategy& ParsedTools::GridRefinement::get_strategy ( ) const
inline

Get the current strategy object.

Definition at line 167 of file grid_refinement.h.

References strategy.

◆ get_n_refinement_cycles()

const unsigned int& ParsedTools::GridRefinement::get_n_refinement_cycles ( ) const
inline

Get the total number of refinemt cycles.

Definition at line 176 of file grid_refinement.h.

References n_refinement_cycles.

◆ get_refinement_cycles()

std_cxx20::ranges::iota_view<unsigned int, unsigned int> ParsedTools::GridRefinement::get_refinement_cycles ( ) const
inline

Get Return an object that can be thought of as an array containing all indices from zero to n_refinement_cycles.

Definition at line 186 of file grid_refinement.h.

References n_refinement_cycles.

◆ estimate_mark_refine() [1/2]

template<int dim, int spacedim, typename VectorType , typename Tria >
void ParsedTools::GridRefinement::estimate_mark_refine ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof_handler,
const VectorType &  solution,
Tria &  tria 
) const

Perform all the steps of the ESTIMATE-MARK-REFINE cycle.

◆ estimate_mark_refine() [2/2]

template<int dim, int spacedim, typename VectorType , typename Tria >
void ParsedTools::GridRefinement::estimate_mark_refine ( const DoFHandler< dim, spacedim > &  dof_handler,
const VectorType &  solution,
Tria &  tria 
) const

Perform all the steps of the ESTIMATE-MARK-REFINE cycle.

◆ limit_levels()

template<int dim, int spacedim>
void ParsedTools::GridRefinement::limit_levels ( Triangulation< dim, spacedim > &  tria) const
private

Make sure that the refinement level is kept between min_level and max_level.

Member Data Documentation

◆ n_refinement_cycles

unsigned int ParsedTools::GridRefinement::n_refinement_cycles
private

Definition at line 220 of file grid_refinement.h.

Referenced by get_n_refinement_cycles(), and get_refinement_cycles().

◆ strategy

RefinementStrategy ParsedTools::GridRefinement::strategy
private

Definition at line 221 of file grid_refinement.h.

Referenced by get_strategy().

◆ top_parameter

double ParsedTools::GridRefinement::top_parameter
private

Definition at line 222 of file grid_refinement.h.

◆ bottom_parameter

double ParsedTools::GridRefinement::bottom_parameter
private

Definition at line 223 of file grid_refinement.h.

◆ max_cells

unsigned int ParsedTools::GridRefinement::max_cells
private

Definition at line 224 of file grid_refinement.h.

◆ min_level

int ParsedTools::GridRefinement::min_level
private

Definition at line 225 of file grid_refinement.h.

◆ max_level

int ParsedTools::GridRefinement::max_level
private

Definition at line 226 of file grid_refinement.h.

◆ estimator_type

std::string ParsedTools::GridRefinement::estimator_type
private

Definition at line 228 of file grid_refinement.h.

◆ optional_estimators

std::map<std::string, std::function<void(Vector<float> &)> > ParsedTools::GridRefinement::optional_estimators
private

Definition at line 230 of file grid_refinement.h.

Referenced by GridRefinement().

◆ component_mask

ComponentMask ParsedTools::GridRefinement::component_mask
private

Definition at line 231 of file grid_refinement.h.

◆ subdomain_id

types::subdomain_id ParsedTools::GridRefinement::subdomain_id
private
Initial value:
=
const types::material_id invalid_material_id

Definition at line 232 of file grid_refinement.h.

◆ material_id

types::material_id ParsedTools::GridRefinement::material_id
private
Initial value:

Definition at line 234 of file grid_refinement.h.


The documentation for this class was generated from the following files: