|
| ElementAnisotropic () |
| Constructor, sets default element definition. More...
|
|
real & | scale () override |
| Reference for element size. More...
|
|
void | set_scale (const real val) override |
| Set the scale for the element.
|
|
real | get_scale () const override |
| Get the scale for the element.
|
|
void | set_anisotropic_ratio (const std::array< real, dim > &ratio) override |
| Set the anisotropic ratio for each reference axis. More...
|
|
std::array< real, dim > | get_anisotropic_ratio () override |
| Get the anisotropic ratio of each reference axis as an array. More...
|
|
real | get_anisotropic_ratio (const unsigned int j) override |
| Get the anisotropic ratio corresponding to the \(j^{th}\) reference axis. More...
|
|
void | set_unit_axis (const std::array< dealii::Tensor< 1, dim, real >, dim > &unit_axis) override |
| setting the (unit) axis direction
|
|
std::array< dealii::Tensor< 1, dim, real >, dim > | get_unit_axis () override |
| getting the (unit) axis direction (array)
|
|
dealii::Tensor< 1, dim, real > | get_unit_axis (const unsigned int j) override |
| getting the (unit) axis direction
|
|
void | set_axis (const std::array< dealii::Tensor< 1, dim, real >, dim > &axis) override |
| setting frame axis j (scaled) at index
|
|
std::array< dealii::Tensor< 1, dim, real >, dim > | get_axis () override |
| getting frame axis j (scaled) at index (array)
|
|
dealii::Tensor< 1, dim, real > | get_axis (const unsigned int j) override |
| getting frame axis j (scaled) at index
|
|
dealii::Tensor< 2, dim, real > | get_metric () override |
| get metric value at index
|
|
dealii::Tensor< 2, dim, real > | get_inverse_metric () override |
| get inverse metric value
|
|
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) override |
| Set anisotropy from reconstructed directional derivatives. More...
|
|
void | apply_anisotropic_limit (const real anisotropic_ratio_min, const real anisotropic_ratio_max) |
| Limits the anisotropic ratios to a given bandwidth. More...
|
|
virtual | ~Element ()=default |
| Destructor.
|
|
template<typename DoFCellAccessorType > |
void | set_cell (const DoFCellAccessorType &cell) |
| Set the Element based on the input cell (from current mesh) More...
|
|
template<int dim, typename real>
class PHiLiP::GridRefinement::ElementAnisotropic< dim, real >
Anisotropic element class.
Specialization of the element type for the anisotropic remeshing case. Stores decomposed frame field axes (size, orientation and anisotropy). A collection of elements is contained in the corresponding Field class.
Definition at line 287 of file field.h.
template<int dim, typename real >
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.
Implements PHiLiP::GridRefinement::Element< dim, real >.
Definition at line 533 of file field.cpp.
template<int dim, typename real >
void PHiLiP::GridRefinement::ElementAnisotropic< dim, real >::set_anisotropy |
( |
const std::array< real, dim > & |
derivative_value, |
|
|
const std::array< dealii::Tensor< 1, dim, real >, dim > & |
derivative_direction, |
|
|
const unsigned int |
order |
|
) |
| |
|
overridevirtual |