|
T | limit (const T &phi_upwind, const T &, const VectorValue< T > *grad_phi_upwind, const VectorValue< T > *, const RealVectorValue &, const Real &max_value, const Real &min_value, const FaceInfo *fi, const bool &fi_elem_is_upwind) const override final |
| This method overrides the pure virtual limit method in the base Limiter class. More...
|
|
bool | constant () const override final |
|
InterpMethod | interpMethod () const override final |
|
| VenkatakrishnanLimiter ()=default |
|
virtual T | limit (const T &phi_upwind, const T &phi_downwind, const libMesh::VectorValue< T > *grad_phi_upwind, const libMesh::VectorValue< T > *grad_phi_downwind, const RealVectorValue &dCD, const Real &max_value, const Real &min_value, const FaceInfo *fi, const bool &fi_elem_is_upwind) const =0 |
| This method computes the flux limiting ratio based on the provided scalar values, gradient vectors, direction vector, maximum and minimum allowable values, and geometric information from the face and cell centroids. More...
|
|
T | operator() (const T &phi_upwind, const T &phi_downwind, const libMesh::VectorValue< T > *grad_phi_upwind, const RealVectorValue &dCD) const |
|
T | operator() (const T &phi_upwind, const T &phi_downwind, const VectorValue< T > *grad_phi_upwind, const VectorValue< T > *grad_phi_downwind, const RealVectorValue &dCD, const Real &max_value, const Real &min_value, const FaceInfo *fi, const bool &fi_elem_is_upwind) const |
|
T | rf_grad (const VectorValue< T > *grad_phi_upwind, const VectorValue< T > *grad_phi_downwind, const RealVectorValue &dCD) const |
|
T | rf_minmax (const T &phi_upwind, const VectorValue< T > *grad_phi_upwind, const Real &max_value, const Real &min_value, const FaceInfo *fi, const bool &fi_elem_is_upwind) const |
|
template<typename T>
class Moose::FV::VenkatakrishnanLimiter< T >
The Venkatakrishnan limiter is derived from the following equations:
- Calculation of the face delta:
\[ \Delta_{\text{face}} = \nabla \phi_{\text{upwind}} \cdot (\mathbf{x}_{\text{face}} - \mathbf{x}_{\text{cell}}) \]
- Calculation of the deltas for maximum and minimum values relative to the upwind value with a small perturbation to avoid division by zero:
\[ \Delta_{\text{max}} = \phi_{\text{max}} - \phi_{\text{upwind}} + 1e-10 \]
\[ \Delta_{\text{min}} = \phi_{\text{min}} - \phi_{\text{upwind}} + 1e-10 \]
- Calculation of the gradient ratio coefficient ( r_f ):
\[ r_f = \begin{cases} \frac{\Delta_{\text{face}}}{\Delta_{\text{max}}} & \text{if } \Delta_{\text{face}} \geq 0 \\ \frac{\Delta_{\text{face}}}{\Delta_{\text{min}}} & \text{otherwise} \end{cases} \]
- Venkatakrishnan limiter formula (Venkatakrishnan, 1993):
\[ \beta(r_f) = \frac{2 r_f + 1.0}{r_f (2 r_f + 1.0) + 1.0} \]
- Template Parameters
-
T | The data type of the scalar values and the return type. |
Definition at line 48 of file VenkatakrishnanLimiter.h.
This method computes the flux limiting ratio based on the provided scalar values, gradient vectors, direction vector, maximum and minimum allowable values, and geometric information from the face and cell centroids.
It must be overridden by any derived class implementing a specific limiting strategy.
- Template Parameters
-
T | The data type of the field values and the return type. |
- Parameters
-
phi_upwind | The field value at the upwind location. |
phi_downwind | The field value at the downwind location. |
grad_phi_upwind | Pointer to the gradient of the field at the upwind location. |
grad_phi_downwind | Pointer to the gradient of the field at the downwind location. |
dCD | A constant direction vector representing the direction of the cell face. |
max_value | The maximum allowable value. |
min_value | The minimum allowable value. |
fi | FaceInfo object containing geometric details such as face centroid and cell centroids. |
fi_elem_is_upwind | Boolean indicating if the face info element is upwind. |
- Returns
- The computed flux limiting ratio.
This pure virtual function is intended to be defined in derived classes, which will implement the specific limiting algorithm. Derived classes will provide the actual logic for computing the flux limiting ratio, ensuring that the result adheres to the constraints and properties required by the specific limiting method.
Referenced by Moose::FV::Limiter< T >::operator()().