[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
#include <field.h>
Public Types | |
using | DoFHandlerType = dealii::DoFHandler< dim > |
Associated DofHandler type. | |
Public Member Functions | |
virtual | ~Field ()=default |
Destructor. | |
virtual void | reinit (const unsigned int size)=0 |
reinitialize the internal data structure | |
virtual unsigned int | size () const =0 |
returns the internal vector size | |
virtual real & | scale (const unsigned int index)=0 |
reference for element size | |
virtual void | set_scale (const unsigned int index, const real val)=0 |
Sets scale for the specified element index. | |
virtual real | get_scale (const unsigned int index) const =0 |
Gets scale for the specified element index. | |
void | set_scale_vector_dealii (const dealii::Vector< real > &vec) |
Sets scale for all elements from a scale vector (dealii::Vector) | |
void | set_scale_vector (const std::vector< real > &vec) |
Sets scale for all elements from a scale vector(std::vector) | |
dealii::Vector< real > | get_scale_vector_dealii () const |
Gets scale for all elements from a scale vector (dealii::Vector) | |
std::vector< real > | get_scale_vector () const |
Gets scale fpr all elements from a scale vector (std::vector) | |
virtual void | set_anisotropic_ratio (const unsigned int index, const std::array< real, dim > &ratio)=0 |
Sets anisotropic ratio for the specified element index (if anisotropic) | |
virtual std::array< real, dim > | get_anisotropic_ratio (const unsigned int index)=0 |
Gets the anisotropic ratio for the specified element index. | |
virtual real | get_anisotropic_ratio (const unsigned int index, const unsigned int j)=0 |
Gets the anisotropic ratio for the specified element index along axis j. | |
dealii::Vector< real > | get_max_anisotropic_ratio_vector_dealii () |
Gets the vector of largest anisotropic ratio for each element (dealii::Vector) | |
virtual void | set_unit_axis (const unsigned int index, const std::array< dealii::Tensor< 1, dim, real >, dim > &unit_axis)=0 |
Sets the group of dim (unit) axis directions for the specified element index. | |
virtual std::array< dealii::Tensor< 1, dim, real >, dim > | get_unit_axis (const unsigned int index)=0 |
Gets the group of dim (unit) axis direction for the specified element index. | |
virtual dealii::Tensor< 1, dim, real > | get_unit_axis (const unsigned int index, const unsigned int j)=0 |
Sets the j^th (unit) axis direction for the specified element index. | |
virtual void | set_axis (const unsigned int index, const std::array< dealii::Tensor< 1, dim, real >, dim > &axis)=0 |
Sets the j^th (scale) axis direction for the specified element index. | |
virtual std::array< dealii::Tensor< 1, dim, real >, dim > | get_axis (const unsigned int index)=0 |
Gets the frame (dim scaled axis vectors) for the specified element index. | |
virtual dealii::Tensor< 1, dim, real > | get_axis (const unsigned int index, const unsigned int j)=0 |
Gets the j^th frame component (scaled axis vector) for the specified element index. | |
void | set_axis_vector (const std::vector< std::array< dealii::Tensor< 1, dim, real >, dim >> &vec) |
Sets the frame (dim scaled axis vectors) for each element (std::vector) | |
std::vector< std::array< dealii::Tensor< 1, dim, real >, dim > > | get_axis_vector () |
Gets the frame (dim scaled axis vectors) for each element (std::vector) | |
std::vector< dealii::Tensor< 1, dim, real > > | get_axis_vector (const unsigned int j) |
Sets the j^th frame component (scaled axis vector) for each element (std::vector) | |
virtual dealii::Tensor< 2, dim, real > | get_metric (const unsigned int index)=0 |
Gets the anisotropic metric tensor, \(M\), for the specified element index. | |
std::vector< dealii::Tensor< 2, dim, real > > | get_metric_vector () |
Gets the anisotropic metric tensor, \(M\), for each element (std::vector) | |
virtual dealii::Tensor< 2, dim, real > | get_inverse_metric (const unsigned int index)=0 |
Gets the inverse metric tensor, \(M^{-1}\), for the specified element index. | |
std::vector< dealii::Tensor< 2, dim, real > > | get_inverse_metric_vector () |
Gets the inverse metric tensor, \(M^{-1}\), for each element (std::vector) | |
dealii::SymmetricTensor< 2, dim, real > | get_quadratic_metric (const unsigned int index) |
Gets the quadratic Riemannian metric, \(\mathcal{M} = M^T M\), for the specified element index. | |
std::vector< dealii::SymmetricTensor< 2, dim, real > > | get_quadratic_metric_vector () |
Gets the quadratic Riemannian metric, \(\mathcal{M} = M^T M\), for each element (std::vector) | |
dealii::SymmetricTensor< 2, dim, real > | get_inverse_quadratic_metric (const unsigned int index) |
Gets the inverse quadratic Riemannian metric used with BAMG, \(\mathcal{M}^{-1}\), for a specified element index. | |
std::vector< dealii::SymmetricTensor< 2, dim, real > > | get_inverse_quadratic_metric_vector () |
Gets the inverse quadratic Riemannian metric used with BAMG, \(\mathcal{M}^{-1}\), for each element (std::vector) | |
virtual void | set_cell (const DoFHandlerType &dof_handler)=0 |
Assigns the existing field based on an input DoFHandlerType. | |
virtual void | set_anisotropy (const dealii::DoFHandler< dim > &dof_handler, const std::vector< std::array< real, dim >> &derivative_value, const std::vector< std::array< dealii::Tensor< 1, dim, real >, dim >> &derivative_direction, const int relative_order)=0 |
Compute anisotropic ratio from directional derivatives. More... | |
virtual void | apply_anisotropic_limit (const real anisotropic_ratio_min, const real anisotropic_ratio_max)=0 |
Globally limit anisotropic ratio to a specified range. | |
virtual std::ostream & | serialize (std::ostream &os) const =0 |
Performs internal call for writing the field description to an ostream. | |
Friends | |
std::ostream & | operator<< (std::ostream &os, const Field< dim, real > &field) |
Performs output to ostream using internal serialization. | |
Field element class.
This class describes the anisotropic size specification for continuous mesh remeshing methods. Each DG mesh element is associated with an isotropic or anisotropic description for the ideal target element at each point during remeshing. Provides wraper to access these values to output to mesh generators without knowledge of mesh type. Each local element is stored in ElementIsotropic or ELementAnisotropic specified above and the FieldInternal class is used to differentiate.
|
pure virtual |
Compute anisotropic ratio from directional derivatives.
Uses Dolejsi's anisotropy method based on reconstructed \(p+1\) directional derivatives. Derivatives are obtained in reconstruct_poly.cpp.
Implemented in PHiLiP::GridRefinement::FieldInternal< dim, real, ElementType >.