|
[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Internal Field element class. More...
#include <field.h>
Public Member Functions | |
| void | reinit (const unsigned int size) |
| reinitialize the internal data structure | |
| unsigned int | size () const |
| returns the internal vector size | |
| real & | scale (const unsigned int index) override |
| reference for element size of specified element index | |
| void | set_scale (const unsigned int index, const real val) override |
| Sets scale for the specified element index. | |
| real | get_scale (const unsigned int index) const override |
| Gets scale for the specified element index. | |
| void | set_anisotropic_ratio (const unsigned int index, const std::array< real, dim > &ratio) override |
| Sets anisotropic ratio for the specified element index (if anisotropic) | |
| std::array< real, dim > | get_anisotropic_ratio (const unsigned int index) override |
| Gets the anisotropic ratio for the specified element index. | |
| real | get_anisotropic_ratio (const unsigned int index, const unsigned int j) override |
| Gets the anisotropic ratio for the specified element index along axis j. | |
| void | set_unit_axis (const unsigned int index, const std::array< dealii::Tensor< 1, dim, real >, dim > &unit_axis) override |
| Sets the group of dim (unit) axis directions for the specified element index. | |
| std::array< dealii::Tensor< 1, dim, real >, dim > | get_unit_axis (const unsigned int index) override |
| Sets the group of dim (unit) axis directions for the specified element index. | |
| dealii::Tensor< 1, dim, real > | get_unit_axis (const unsigned int index, const unsigned int j) override |
| Gets the group of dim (unit) axis direction for the specified element index. | |
| void | set_axis (const unsigned int index, const std::array< dealii::Tensor< 1, dim, real >, dim > &axis) override |
| Sets the j^th (unit) axis direction for the specified element index. | |
| std::array< dealii::Tensor< 1, dim, real >, dim > | get_axis (const unsigned int index) override |
| Gets the frame (dim scaled axis vectors) for the specified element index. | |
| dealii::Tensor< 1, dim, real > | get_axis (const unsigned int index, const unsigned int j) override |
| Gets the j^th frame component (scaled axis vector) for the specified element index. | |
| dealii::Tensor< 2, dim, real > | get_metric (const unsigned int index) override |
| Gets the anisotropic metric tensor, \(M\), for the specified element index. | |
| dealii::Tensor< 2, dim, real > | get_inverse_metric (const unsigned int index) override |
| Gets the inverse metric tensor, \(M^{-1}\), for the specified element index. | |
| void | set_cell (const typename Field< dim, real >::DoFHandlerType &dof_handler) override |
| 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) override |
| Compute anisotropic ratio from directional derivatives. More... | |
| void | apply_anisotropic_limit (const real anisotropic_ratio_min, const real anisotropic_ratio_max) override |
| Globally limit anisotropic ratio to a specified range. | |
| std::ostream & | serialize (std::ostream &os) const override |
| Performs internal call for writing the field description to an ostream. | |
Public Member Functions inherited from PHiLiP::GridRefinement::Field< dim, real > | |
| virtual | ~Field ()=default |
| Destructor. | |
| 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) | |
| dealii::Vector< real > | get_max_anisotropic_ratio_vector_dealii () |
| Gets the vector of largest anisotropic ratio for each element (dealii::Vector) | |
| 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) | |
| std::vector< dealii::Tensor< 2, dim, real > > | get_metric_vector () |
| Gets the anisotropic metric tensor, \(M\), for each element (std::vector) | |
| 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. | |
Private Attributes | |
| std::vector< ElementType > | field |
| Internal vector storage of element data for each index. | |
Additional Inherited Members | |
Public Types inherited from PHiLiP::GridRefinement::Field< dim, real > | |
| using | DoFHandlerType = dealii::DoFHandler< dim > |
| Associated DofHandler type. | |
Internal Field element class.
This class reimplements the above generalized function with internal handling of the ElementType. This prevents i/o and certain adaption functions from needing to directly check the type of field which is being used. Internally stores the target output mesh description from the continuous space as a discrete vector associated with each element.
|
overridevirtual |
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.
Implements PHiLiP::GridRefinement::Field< dim, real >.