[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Fixed-Fraction Grid Refinement Class. More...
#include <grid_refinement_fixed_fraction.h>
Public Member Functions | |
void | refine_grid () override |
Perform call to the grid refinement object of choice. More... | |
![]() | |
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 () |
Based on error indicator, performs fixed-fraction flagging for element splitting. More... | |
void | refine_grid_p () |
Based on error indicator, performs fixed-fraction flagging for polynomial order enrichment. More... | |
void | refine_grid_hp () |
(NOT IMPLEMENTED) Based on error and smoothness indicator, perform fixed-fraction flagging and decision between element splitting and polynomial order enrichment More... | |
void | refine_boundary_h () |
Flags the domain boundaries for refinement. More... | |
void | smoothness_indicator () |
(NOT IMPLEMENTED) Computes smoothness indicator for \(hp\)-refinement decisions More... | |
void | anisotropic_h () |
Performs anisotropic adaptation of the mesh. More... | |
void | anisotropic_h_jump_based () |
Sets anisotropic refinement flags for cells based on discontinuity between neighbors. More... | |
void | anisotropic_h_reconstruction_based () |
Sets anisotropic refinement flags for cells based on directional derivatives reconstructed along the chord lines. More... | |
void | error_indicator () |
Performs call to proper error indcator function based on type of error_indicator parameter. More... | |
void | error_indicator_error () |
Compute error indicator based on Lq norm relative to exact manufactured solution. More... | |
void | error_indicator_hessian () |
Compute error indicator based on reconstructed \(p+1\) directional derivatives. More... | |
void | error_indicator_residual () |
(NOT IMPLEMENTED) Compute error indicator based on fine grid residual distribution More... | |
void | error_indicator_adjoint () |
Compute error indicator for goal-oriented approach using a dual-weighted residual. More... | |
std::vector< std::pair< dealii::Vector< real >, std::string > > | output_results_vtk_method () override |
Output refinement method dependent results. More... | |
![]() | |
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. | |
Protected Attributes | |
dealii::Vector< real > | indicator |
dealii::Vector< real > | smoothness |
![]() | |
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... | |
Additional Inherited Members | |
![]() | |
using | ErrorIndicatorEnum = PHiLiP::Parameters::GridRefinementParam::ErrorIndicator |
Fixed-Fraction Grid Refinement Class.
This class offers methods to perform fixed-fraction style refinement using built in Deal.II refinement methods. These methods work by first flagging a chosen fraction of the mesh with largest error indicator to be refined. In \(h\)-refinement, this is applied by splitting each cell into a set of subcells (which may potentially be anisotropic). Methods for \(p\)-refinement are also supported which locally increment the polynomial order of the discretization. \(hp\)-refinement is not fully supported, but some placeholder functions have been included to facilitate this in the future. Note: Some functionality and behaviour will vary depending on the choice of triangulation used for both the refinement step and in other parts of the code.
Definition at line 27 of file grid_refinement_fixed_fraction.h.
|
protected |
Performs anisotropic adaptation of the mesh.
Based on the choice of aniso_indicator, will call either jump_based or reconstruction_based anisotropy methods. Note: anisotropic adaptation is not currently applicable to all deal.II mesh types.
Definition at line 173 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Sets anisotropic refinement flags for cells based on discontinuity between neighbors.
If the amount of jump across pairs of face in a given element axes are sufficiently above the average of all faces for the cell, then that axis will be flagged for anisotropic refinement instead of regular refinement. The anisotropic_threshold_ratio parameter controls how strong an anisotropic behavior triggers the response. This option is only considered on cells flagged isotropically such that the refinement strategy is modified before the application of the splitting. This is based on Step 30 of the Deal.II examples: https://www.dealii.org/current/doxygen/deal.II/step_30.html
Definition at line 188 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Sets anisotropic refinement flags for cells based on directional derivatives reconstructed along the chord lines.
Uses a polynomial reconstruction of the high-order \(p+1\) derivatives to estimate the effectivness of refinement along a single cell axis instead of isotropically. The indicator along each cell chord (the vector from opposing face center to opposing face center) and compared relative to the anisotropic_threshold_ratio needed to alter the flagging of the cell to be split along a single axis instead. Unlike continuous refinement methods using reconstruction the effectiveness of anisotropy is limited by the initial cell orientation and dimensions.
Definition at line 318 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Performs call to proper error indcator function based on type of error_indicator parameter.
Each of these functions operates by filling the indicator member variable vector with a value for each cell in the mesh. These inidcators are then used to decide which cells of the mesh will be flagged for coarsening and refinement respectively based on the established fractions.
Definition at line 388 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Compute error indicator for goal-oriented approach using a dual-weighted residual.
Uses the sensitivity of the local flow solution relative to the functional of interest obtained through the solution of the adjoint problem to weight the fine grid residual. Together this forms the "dual-weighted residual (DWR)" and provides an estimate of the contribution of local errors in the flow to the global goal-oriented functional evaluation erorrs. This can be written globally as:
\[ |\mathcal{J}(u) - \mathcal{J}_h(u_h^H)|\approx-\phi_h^T\mathcal{R}_h(u_h^H) \]
where locally, the DWR is evaluated by a sum over fine grid elements:
\[ E_i = |\sum_{\Omega_k \in \Omega_i} {(\psi_h)^T_k (\mathcal{R}_h(u_h^H))_k}| \]
which has been implemented as part of the adjoint class using a fine grid adjoint solution.
Definition at line 557 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Compute error indicator based on Lq norm relative to exact manufactured solution.
For debugging and testing, uses the manufactured solution function to determine innacuracies in the current results. Not useful for any actual refinement as it requires knowledge about the exact problem solution, \(u(\boldsymbol{x}\):
\[ E_i = \int_{\Omega_i} {|u_h(\boldsymbol{x}) - u(\boldsymbol{x})|^{Lq} \mathrm{d}\boldsymbol{x}} \]
Evaluated for each element \(\Omega_i\) in the mesh.
Definition at line 406 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Compute error indicator based on reconstructed \(p+1\) directional derivatives.
Uses the reconstructed enriched solution space to determine which areas of the mesh result in the largest discretization error. It is similar in strategy to existing hessian-based and similar high-order approximations. However, in the current scope it is applied to flag cells for splitting. This is done based on the fact the error approximates the \(p+1\) solution:
\[ E_i = u_{h,p+1}(\boldsymbol{x}+h\boldsymbol{\xi}) - u_{h,p}(\boldsymbol{x}+h\boldsymbol{\xi}) = D_{\boldsymbol{\xi}}^{p+1} u(x) h^{p+1} \]
based on the taylor series expansion where \(\boldsymbol{\xi}\) is the direction vector used in the largest directional derivative
Definition at line 457 of file grid_refinement_fixed_fraction.cpp.
|
protected |
(NOT IMPLEMENTED) Compute error indicator based on fine grid residual distribution
By reconstructing the current solution on a finer grid and evaluating the magnitude of the local residuals, this strategy approximates how innacurate the current approximation would be in a smoother space or after the local mesh is refined. It is similar to a goal-oriented approach without adjoint weighting:
\[ E_i = |\mathcal{R}_h(u_h^H)| \]
where \(u_h^H\) is the prolongation operator to the fine mesh.
Definition at line 507 of file grid_refinement_fixed_fraction.cpp.
|
overrideprotectedvirtual |
Output refinement method dependent results.
For the current class this includes the error indicator and smoothness indicator distribtuions.
Implements PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >.
Definition at line 572 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Flags the domain boundaries for refinement.
Used for testing the effect this resolution has on the measurement of the functional output value for boundary integral based functionals.
Definition at line 144 of file grid_refinement_fixed_fraction.cpp.
|
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.
For fixed-fraction refinement, this function first computes the error-indicator and optionally smoothness indicator (if used for anisotropy or \(hp\)-refinement) then calls proper refinement function to flag and modify the grid in the desired manner (locally splitting cells or change polynomial orderes). Also automatically transfers the solution onto the new embedded mesh.
See subclass functions for details of refinement types.
Implements PHiLiP::GridRefinement::GridRefinementBase< dim, nstate, real, MeshType >.
Definition at line 15 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Based on error indicator, performs fixed-fraction flagging for element splitting.
Based on the distribution of the error indicator, the fraction of cells with the largest values will be marked for refinement and the fraction with the lowest values will be marked for coarsening. At the end of the refine_grid function, when the refinement is executed refined cells will be split into \(2^{dim}\) subcells along with any neighbors needed to maintain 2:1 face connectivity. Coarsened cells will be recombined if neighbor cells have also been marked in such a way that they can be merged. Usually coarsening is less consistently applied and it may be simpler to only target refinements.
Definition at line 87 of file grid_refinement_fixed_fraction.cpp.
|
protected |
(NOT IMPLEMENTED) Based on error and smoothness indicator, perform fixed-fraction flagging and decision between element splitting and polynomial order enrichment
After first flagging the cells using the refine_grid_h function, this method loops back over the grid and selectively changes the local cell to target \(p\)-refinement if the smoothness indicator is above a specified tolerance. This offers better overall convergence in smooth areas of the flow, but, due to Gibb's phenomenon, \(h\)-refinement is better suited for areas with discontinuities (such as shocks).
Definition at line 120 of file grid_refinement_fixed_fraction.cpp.
|
protected |
Based on error indicator, performs fixed-fraction flagging for polynomial order enrichment.
First uses the refine_grid_h function to perform flagging of worst portion of cells then loops over the mesh and converts each of these flags to instead increment the target finite element polynomial order for the local cell instead. No coarsening is currently supported.
Definition at line 106 of file grid_refinement_fixed_fraction.cpp.
|
protected |
(NOT IMPLEMENTED) Computes smoothness indicator for \(hp\)-refinement decisions
Due to Gibbs phenomenon, even with a high-order method convergence degrades in regions with discontinuities. Therefore, for capturing shocks and other phenomenon \(h\)-refinement is more suitable. However, in smoother areas of the problem, \(p\)-refinement is capable of capturing large areas with fewer elements offering potential spectral convergence. Therefore, this placeholder function is intended to approximate the smoothness of the flow for choosing between the refinement types.
Definition at line 161 of file grid_refinement_fixed_fraction.cpp.