[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Base Grid Refinement Class. More...
#include <grid_refinement.h>
Public Member Functions | |
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. | |
virtual void | refine_grid ()=0 |
Perform call to the grid refinement object of choice. More... | |
void | output_results_vtk (const unsigned int iref) |
Write information about the grid refinement step to a .vtk file. More... | |
Protected Types | |
using | ErrorIndicatorEnum = PHiLiP::Parameters::GridRefinementParam::ErrorIndicator |
Protected Member Functions | |
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. | |
virtual std::vector< std::pair< dealii::Vector< real >, std::string > > | output_results_vtk_method ()=0 |
Output refinement method dependent results. More... | |
Protected Attributes | |
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... | |
Base Grid Refinement Class.
This class provides access to the basic refinement control methods implemented for uniform, fixed-fraction and continuous style methods in the associated *.cpp files. Although this class contains no refinement functionality of its own, the virtual functions implemented from here provide a uniform interface for adapting the grid based on a variety of h-, p- and hp- style refinement techniques with indicators from the exact error (manufactured solution), feature-based (generalization of hessian based methods for high order), the local residual distribution on the fine grid and goal-oriented adjoint-based techniques.
Additionally, this class contains functionality for writing a description of the current grid refinement object to a .vtk file with additional refinement information passed from the subclass implementations.
See the related parameter object PHiLiP::Parameters::GridRefinementParam for more information about the various options and controls availible.
Note: This class templated on the mesh type as anisotropic fixed-fraction splitting is not availible in parralel at this time.
Definition at line 42 of file grid_refinement.h.
void PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >::output_results_vtk | ( | const unsigned int | iref | ) |
Write information about the grid refinement step to a .vtk file.
Includes various information about the mesh, solution, error indicators, target refinements (both h- and p-), functional solution, physics, etc.
also provides interface for subclasses to output additional visualization fields.
Definition at line 40 of file grid_refinement.cpp.
|
protectedpure virtual |
Output refinement method dependent results.
This class is overridden in the subclasses with any additional visualization fields.
Implemented in PHiLiP::GridRefinement::GridRefinement_FixedFraction< dim, nstate, real, MeshType >, PHiLiP::GridRefinement::GridRefinement_Continuous< dim, nstate, real, MeshType >, and PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >.
|
pure virtual |
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.
Implemented in PHiLiP::GridRefinement::GridRefinement_Continuous< dim, nstate, real, MeshType >, PHiLiP::GridRefinement::GridRefinement_FixedFraction< dim, nstate, real, MeshType >, and PHiLiP::GridRefinement::GridRefinement_Uniform< dim, nstate, real, MeshType >.
|
protected |
update flags needed at face points
Definition at line 186 of file grid_refinement.h.
|
protected |
update flags needed at neighbor's face points
Definition at line 195 of file grid_refinement.h.
|
protected |
Triangulation object of templated mesh type.
Note: anisotropic, p- type and other refinements may not work in all cases.
Definition at line 169 of file grid_refinement.h.
|
protected |
update flags needed at volume points
Definition at line 178 of file grid_refinement.h.