[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
#include <anisotropic_mesh_adaptation.h>
Public Member Functions | |
AnisotropicMeshAdaptation (std::shared_ptr< DGBase< dim, real, MeshType > > dg_input, const real _norm_Lp, const real _complexity, const bool _use_goal_oriented_approach=false) | |
Constructor. | |
void | adapt_mesh () |
Function which adapts mesh and loads in new mesh. | |
Private Types | |
using | VectorType = dealii::LinearAlgebra::distributed::Vector< double > |
Alias for dealii's parallel distributed vector. | |
Private Member Functions | |
dealii::Tensor< 2, dim, real > | get_positive_definite_tensor (const dealii::Tensor< 2, dim, real > &input_tensor) const |
Returns positive tensor from an input tensor by taking absolute of eigenvalues. | |
void | compute_cellwise_optimal_metric () |
Computes optimal metric depending on goal oriented or feature based approach. | |
void | compute_abs_hessian () |
Computes hessian using the input coefficients, which can be a solution sensor or (for goal oriented approach) convective flux. More... | |
void | initialize_cellwise_metric_and_hessians () |
Initializes cellwise metric and hessian to zero tensors. | |
void | compute_feature_based_hessian () |
Computes feature based hessian (i.e. hessian of the solution). | |
void | compute_goal_oriented_hessian () |
Computes pseudo Hessian for the goal oriented approach. | |
void | change_p_degree_and_interpolate_solution (const unsigned int poly_degree) |
Change the polynomial order and interpolate solution. | |
void | reconstruct_p2_solution () |
Reconstructs p2 solution after interpolation. More... | |
unsigned int | get_iquad_near_cellcenter (const dealii::Quadrature< dim > &volume_quadrature) const |
Returns quadrature number of a point which is closest to the reference cell's center. | |
void | get_flux_at_support_pts (std::vector< std::array< dealii::Tensor< 1, dim, real >, nstate >> &flux_at_support_pts, const dealii::FEValues< dim, dim > &fe_values_input, const std::vector< dealii::types::global_dof_index > &dof_indices, typename dealii::DoFHandler< dim >::active_cell_iterator cell) const |
Returns flux coeffs by evaluating flux at support points of fe. More... | |
Private Attributes | |
std::shared_ptr< DGBase< dim, real, MeshType > > | dg |
Shared pointer to DGBase. | |
const bool | use_goal_oriented_approach |
Flag to use goal oriented approach. It is set to false by default. | |
std::vector< dealii::Tensor< 2, dim, real > > | cellwise_hessian |
Stores hessian in each cell. | |
const real | normLp |
Lp Norm w.r.t. which the anlytical optimization is done. | |
const real | complexity |
Analogue of number of vertices/elements in continuous mesh framework. | |
MPI_Comm | mpi_communicator |
Alias for MPI_COMM_WORLD. | |
dealii::ConditionalOStream | pcout |
std::cout only on processor #0. | |
int | mpi_rank |
Processor# of current processor. | |
int | n_mpi |
Total no. of processors. | |
const unsigned int | initial_poly_degree |
Stores initial polynomial degree. | |
std::shared_ptr< Physics::PhysicsBase< dim, nstate, real > > | pde_physics_double |
Contains the physics of the PDE with real type. | |
std::shared_ptr< Functional< dim, nstate, real, MeshType > > | functional |
Functional to evaluate the adjoint for goal oriented anisotropic meshing. | |
std::vector< dealii::Tensor< 2, dim, real > > | cellwise_optimal_metric |
Stores optimal metric in each cell. | |
const unsigned int | max_dofs_per_cell |
Stores max dofs per cell to initialize dof_indices. | |
Performs anisotropic mesh adaptation with an optimal metric for P1 solution. Implements the optimal metric field derived from continuous optimization framework. See papers from INRIA for further details: Feature based : Loseille, A. and Alauzet, F. "Continuous mesh framework part I: well-posed continuous interpolation error.", 2011. Goal oriented: Loseille, A., Dervieux, A., and Alauzet, F. "Fully anisotropic goal-oriented mesh adaptation for 3d steady Euler equations.", 2010.
Definition at line 22 of file anisotropic_mesh_adaptation.h.
|
private |
Computes hessian using the input coefficients, which can be a solution sensor or (for goal oriented approach) convective flux.
This function is called by compute_optimal_metric().
Definition at line 198 of file anisotropic_mesh_adaptation.cpp.
|
private |
Returns flux coeffs by evaluating flux at support points of fe.
Referenced by flux[idof][istate][idim]
Definition at line 422 of file anisotropic_mesh_adaptation.cpp.
|
private |
Reconstructs p2 solution after interpolation.
Currenlty done using 1 linear solve of the implicit system.
Definition at line 251 of file anisotropic_mesh_adaptation.cpp.