[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 | VertexList = std::array< dealii::Tensor< 1, dim, real >, dealii::GeometryInfo< dim >::vertices_per_cell > |
Type alias for array of vertex coordinates for element. | |
using | ChordList = std::array< dealii::Tensor< 1, dim, real >, dim > |
Type alias for array of chord veectors (face center to face center) in Deal.II ordering. | |
Public Member Functions | |
virtual | ~Element ()=default |
Destructor. | |
virtual real & | scale ()=0 |
Reference for element size. More... | |
virtual void | set_scale (const real val)=0 |
Set the scale for the element. | |
virtual real | get_scale () const =0 |
Get the scale for the element. | |
virtual void | set_anisotropic_ratio (const std::array< real, dim > &ratio)=0 |
Set the anisotropic ratio for each reference axis. More... | |
virtual std::array< real, dim > | get_anisotropic_ratio ()=0 |
Get the anisotropic ratio of each reference axis as an array. More... | |
virtual real | get_anisotropic_ratio (const unsigned int j)=0 |
Get the anisotropic ratio corresponding to the \(j^{th}\) reference axis. More... | |
virtual void | set_unit_axis (const std::array< dealii::Tensor< 1, dim, real >, dim > &unit_axis)=0 |
Set the unit axes of the element. More... | |
virtual std::array< dealii::Tensor< 1, dim, real >, dim > | get_unit_axis ()=0 |
Get unit axis directions as an array. | |
virtual dealii::Tensor< 1, dim, real > | get_unit_axis (const unsigned int j)=0 |
Get the \(j^{th}\) reference axis. | |
virtual void | set_axis (const std::array< dealii::Tensor< 1, dim, real >, dim > &axis)=0 |
Set the scaled local frame axes based on vector set (length and direction) More... | |
virtual std::array< dealii::Tensor< 1, dim, real >, dim > | get_axis ()=0 |
Get the array of axes of the local frame field \((v_1,v_2,v_3)\). | |
virtual dealii::Tensor< 1, dim, real > | get_axis (const unsigned int j)=0 |
Get the vector corresponding to the \(j^{th}\) frame axis \(v_j\). | |
virtual dealii::Tensor< 2, dim, real > | get_metric ()=0 |
Get metric matrix at point describing mapping from reference element. More... | |
virtual dealii::Tensor< 2, dim, real > | get_inverse_metric ()=0 |
. Get inverse metric matrix for the reference element More... | |
template<typename DoFCellAccessorType > | |
void | set_cell (const DoFCellAccessorType &cell) |
Set the Element based on the input cell (from current mesh) More... | |
virtual void | set_anisotropy (const std::array< real, dim > &derivative_value, const std::array< dealii::Tensor< 1, dim, real >, dim > &derivative_direction, const unsigned int order)=0 |
Set anisotropy from reconstructed directional derivatives. More... | |
virtual void | apply_anisotropic_limit (const real anisotropic_ratio_min, const real anisotropic_ratio_max)=0 |
Limits the anisotropic ratios to a given bandwidth. More... | |
Protected Member Functions | |
ChordList | get_chord_list (const VertexList &vertices) |
Get the chord list from an input set of vertices. | |
virtual void | set_cell_internal (const VertexList &vertices)=0 |
Set element to match geometry of input vertex set. More... | |
Element class.
This provides a base class for incorporating both isotropic (size) and anisotropic (size, anisotropy orientation) in controlling unstructured mesh adaptation methods. Provides functions for setting and accesing the local values for a given cell. A collection of elements is contained in the corresponding Field class (with respective extensions for anisotropy). Note: setting anisotropic properties in the isotropic case will assert(0) and make no changes.
|
pure virtual |
Limits the anisotropic ratios to a given bandwidth.
First finds ratio above max, redistributes length change to maintain constant volume and scale. Then the process is repeated with a lower bound. Note: does nothing in the isotropic case.
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
pure virtual |
Get the anisotropic ratio of each reference axis as an array.
Note: equals 1 for each axis in isotropic case.
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
pure virtual |
Get the anisotropic ratio corresponding to the \(j^{th}\) reference axis.
Note: equals 1 for each axis in isotropic case.
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
pure virtual |
. Get inverse metric matrix for the reference element
In 2D orthogonal case, \(V = 1/h * \mathrm{diag}{\rho,1/\rho} * R(-\theta)\).
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
pure virtual |
Get metric matrix at point describing mapping from reference element.
In 2D orthorgonal case, \(V = [v,w] = h * R(\theta) * \mathrm{diag}{\rho,1/\rho}\). Under transformation, order of axes is maintained with \((i,j,k)\) vectors mapping to \((v_1,v_2,v_3)\)
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
pure virtual |
Reference for element size.
Allows direct read/write of scale of mean element axis. Measure of element (length, area, volume) will be \(scale^{dim}\)
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
pure virtual |
Set the anisotropic ratio for each reference axis.
Requires array of order matching axis definition. Each reference axis will have length \(l = \alpha * scale\). Note: does nothing in the isotropic case.
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
pure virtual |
Set anisotropy from reconstructed directional derivatives.
Based on Dolejsi's method for simplices, uses values obtained from reconstructed polynomial on neighbouring cell patch.
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
pure virtual |
Set the scaled local frame axes based on vector set (length and direction)
Describe the axes of reference parrelogram or parrelopiped at point. Note: does nothing in the isotropic case.
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.
|
inline |
|
protectedpure virtual |
Set element to match geometry of input vertex set.
Vertices describe the tensor product element in Deal.II ordering. Internal function used in handling of set_cell().
|
pure virtual |
Set the unit axes of the element.
Requires an ordered array of dealii::Tensor with vector axes. If none unit values are provided, will be rescaled and factored in to anisotropic ratios. Note: does nothing in the isotropic case.
Implemented in PHiLiP::GridRefinement::ElementAnisotropic< dim, real >, and PHiLiP::GridRefinement::ElementIsotropic< dim, real >.