[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType > Class Template Reference

Uniform Grid Refinement Class. More...

#include <grid_refinement_uniform.h>

Inheritance diagram for PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >:
Collaboration diagram for PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >:

Public Member Functions

void refine_grid () override
 Perform call to the grid refinement object of choice. More...
 
- Public Member Functions inherited from PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >
virtual ~GridRefinementBase ()=default
 Destructor.
 
 GridRefinementBase (PHiLiP::Parameters::GridRefinementParam gr_param_input, std::shared_ptr< PHiLiP::Adjoint< dim, nstate, real, MeshType > > adj_input, std::shared_ptr< PHiLiP::Physics::PhysicsBase< dim, nstate, real > > physics_input)
 Constructor. Stores the adjoint object, physics and parameters.
 
 GridRefinementBase (PHiLiP::Parameters::GridRefinementParam gr_param_input, std::shared_ptr< PHiLiP::DGBase< dim, real, MeshType > > dg_input, std::shared_ptr< PHiLiP::Physics::PhysicsBase< dim, nstate, real > > physics_input, std::shared_ptr< PHiLiP::Functional< dim, nstate, real, MeshType > > functional_input)
 Constructor. Storers the dg object, physics, functional and parameters.
 
 GridRefinementBase (PHiLiP::Parameters::GridRefinementParam gr_param_input, std::shared_ptr< PHiLiP::DGBase< dim, real, MeshType > > dg_input, std::shared_ptr< PHiLiP::Physics::PhysicsBase< dim, nstate, real > > physics_input)
 Constructor. Stores the dg object, physics and parameters.
 
 GridRefinementBase (PHiLiP::Parameters::GridRefinementParam gr_param_input, std::shared_ptr< PHiLiP::DGBase< dim, real, MeshType > > dg_input)
 Constructor. Stores the dg object and parameters.
 
void output_results_vtk (const unsigned int iref)
 Write information about the grid refinement step to a .vtk file. More...
 

Protected Member Functions

void refine_grid_h ()
 Uniform \(h\) grid refinement. More...
 
void refine_grid_p ()
 Uniform \(p\) grid refinement. More...
 
void refine_grid_hp ()
 Uniform \(hp\) grid refinement. More...
 
std::vector< std::pair< dealii::Vector< real >, std::string > > output_results_vtk_method () override
 Output refinement method dependent results. More...
 
- Protected Member Functions inherited from PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >
 GridRefinementBase (PHiLiP::Parameters::GridRefinementParam gr_param_input, std::shared_ptr< PHiLiP::Adjoint< dim, nstate, real, MeshType > > adj_input, std::shared_ptr< PHiLiP::Functional< dim, nstate, real, MeshType > > functional_input, std::shared_ptr< PHiLiP::DGBase< dim, real, MeshType > > dg_input, std::shared_ptr< PHiLiP::Physics::PhysicsBase< dim, nstate, real > > physics_input)
 Delegated constructor which handles the various optional inputs and setup.
 
void output_results_vtk_dg (dealii::DataOut< dim, dealii::DoFHandler< dim >> &data_out, std::shared_ptr< dealii::DataPostprocessor< dim > > &post_processor, dealii::Vector< float > &subdomain, std::vector< unsigned int > &active_fe_indices, dealii::Vector< double > &cell_poly_degree, std::vector< std::string > &residual_names)
 Output refinement results related to the DG object.
 
void output_results_vtk_functional (dealii::DataOut< dim, dealii::DoFHandler< dim >> &data_out)
 Output refinement results related to the functional object.
 
void output_results_vtk_physics (dealii::DataOut< dim, dealii::DoFHandler< dim >> &data_out)
 Output refinement results related to the problem physics.
 
void output_results_vtk_adjoint (dealii::DataOut< dim, dealii::DoFHandler< dim >> &data_out, std::vector< std::string > &dIdw_names_coarse, std::vector< std::string > &adjoint_names_coarse, std::vector< std::string > &dIdw_names_fine, std::vector< std::string > &adjoint_names_fine)
 Output refinement results related to the adjoint object.
 
void output_results_vtk_error (dealii::DataOut< dim, dealii::DoFHandler< dim >> &data_out, dealii::Vector< real > &l2_error_vec)
 Output refinement results related to the solution error.
 

Additional Inherited Members

- Protected Types inherited from PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >
using ErrorIndicatorEnum = PHiLiP::Parameters::GridRefinementParam::ErrorIndicator
 
- Protected Attributes inherited from PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >
PHiLiP::Parameters::GridRefinementParam grid_refinement_param
 Grid refinement parameters.
 
ErrorIndicatorEnum error_indicator_type
 
std::shared_ptr< PHiLiP::Adjoint< dim, nstate, real, MeshType > > adjoint
 Adjoint object (if provided to factory)
 
std::shared_ptr< PHiLiP::Functional< dim, nstate, real, MeshType > > functional
 Functional object (if provided to factory, directly or indirectly)
 
std::shared_ptr< PHiLiP::DGBase< dim, real, MeshType > > dg
 Discontinuous Galerkin object.
 
std::shared_ptr< PHiLiP::Physics::PhysicsBase< dim, nstate, real > > physics
 Problem physics (if provided to factory, directly or indirectly)
 
std::shared_ptr< MeshType > tria
 Triangulation object of templated mesh type. More...
 
unsigned int iteration
 Internal refinement steps iteration counter.
 
MPI_Comm mpi_communicator
 MPI communicator.
 
dealii::ConditionalOStream pcout
 Parallel std::cout that only outputs on mpi_rank==0.
 
const dealii::UpdateFlags volume_update_flags
 update flags needed at volume points More...
 
const dealii::UpdateFlags face_update_flags
 update flags needed at face points More...
 
const dealii::UpdateFlags neighbor_face_update_flags
 update flags needed at neighbor's face points More...
 

Detailed Description

template<int dim, int nstate, typename real, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
class PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >

Uniform Grid Refinement Class.

Simplest form of grid refinement, included primarily as a benchmark for other strategies. This simple class always applies changes in a consistent way throughout the entire mesh, either by subdividing the mesh cells to form smaller cells ($h$-refinement) incrementing the polynomial orders ($p$-refinement) or a combination of the two ($hp$-refinement).

Definition at line 25 of file grid_refinement_uniform.h.

Member Function Documentation

◆ output_results_vtk_method()

template<int dim, int nstate, typename real , typename MeshType >
std::vector< std::pair< dealii::Vector< real >, std::string > > PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >::output_results_vtk_method ( )
overrideprotectedvirtual

Output refinement method dependent results.

No additional results are included for uniform grid refinements.

Implements PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >.

Definition at line 80 of file grid_refinement_uniform.cpp.

◆ refine_grid()

template<int dim, int nstate, typename real , typename MeshType >
void PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >::refine_grid ( )
overridevirtual

Perform call to the grid refinement object of choice.

This will automatically select the proper subclass, error indicator and various refinement types based on the grid refinement parameters passed at setup to the grid refinement factor class.

See subclass functions for details of refinement types.

Implements PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >.

Definition at line 14 of file grid_refinement_uniform.cpp.

◆ refine_grid_h()

template<int dim, int nstate, typename real , typename MeshType >
void PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >::refine_grid_h ( )
protected

Uniform \(h\) grid refinement.

Uniformly refines the mesh with \(h\)-refinement by subdividing all cells.

Definition at line 58 of file grid_refinement_uniform.cpp.

◆ refine_grid_hp()

template<int dim, int nstate, typename real , typename MeshType >
void PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >::refine_grid_hp ( )
protected

Uniform \(hp\) grid refinement.

Performs call to uniform mesh \(h\)-refinement by splitting all cells followed by uniform \(p\)-refinement where the polynomial order of each cell is incremented by 1.

Definition at line 73 of file grid_refinement_uniform.cpp.

◆ refine_grid_p()

template<int dim, int nstate, typename real , typename MeshType >
void PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >::refine_grid_p ( )
protected

Uniform \(p\) grid refinement.

Increments the polynomial order of each cell by 1.

Definition at line 64 of file grid_refinement_uniform.cpp.


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