Kokkos kernel that computes the Zalesak nodal correction factors \(R_i^+\) and \(R_i^-\) from the low-order predictor \(T^L\) and the antidiffusive fluxes. More...
#include <fct_advection_diffusion.hpp>
Public Types | |
| using | GH = fct_detail::GeometryHelper< ScalarT > |
Public Member Functions | |
| void | operator() (const int id, const int x, const int y, const int r) const |
Public Attributes | |
| grid::Grid4DDataScalar< ScalarT > | T_L_ |
| \(T^L\) with ghost layers filled. | |
| grid::Grid5DDataScalar< ScalarT > | antidiff_ |
| \(\tilde{f}_{ij}\): pre-scaled antidiff fluxes. | |
| grid::Grid4DDataScalar< ScalarT > | R_plus_ |
| Output: \(R_i^+\) correction factor. | |
| grid::Grid4DDataScalar< ScalarT > | R_minus_ |
| Output: \(R_i^-\) correction factor. | |
Static Public Attributes | |
| static constexpr int | num_neighbors = GH::num_neighbors |
Kokkos kernel that computes the Zalesak nodal correction factors \(R_i^+\) and \(R_i^-\) from the low-order predictor \(T^L\) and the antidiffusive fluxes.
Per cell \(i\) (summing over all 6 neighbours \(j \in \mathcal{N}(i)\)):
Positive/negative flux sums:
\[ P_i^+ = \sum_{j:\,\tilde{f}_{ij} > 0} \tilde{f}_{ij}, \qquad P_i^- = \sum_{j:\,\tilde{f}_{ij} < 0} \tilde{f}_{ij}. \]
Local extrema of the low-order solution:
\[ T_i^{\max} = \max(T_i^L,\;T_j^L\;\forall j \in \mathcal{N}(i)), \qquad T_i^{\min} = \min(T_i^L,\;T_j^L\;\forall j \in \mathcal{N}(i)). \]
Room to grow/shrink:
\[ Q_i^+ = T_i^{\max} - T_i^L \geq 0, \qquad Q_i^- = T_i^{\min} - T_i^L \leq 0. \]
Correction factors (Zalesak 1979, eq. 13–14):
\[ R_i^+ = \begin{cases} \min\!\left(1,\;\dfrac{Q_i^+}{P_i^+}\right) & P_i^+ > 0 \\ 1 & \text{otherwise} \end{cases}, \qquad R_i^- = \begin{cases} \min\!\left(1,\;\dfrac{Q_i^-}{P_i^-}\right) & P_i^- < 0 \\ 1 & \text{otherwise} \end{cases}. \]
| ScalarT | Floating-point scalar type. |
| using terra::fv::hex::operators::FCTLimiterKernel< ScalarT >::GH = fct_detail::GeometryHelper< ScalarT > |
|
inline |
| grid::Grid5DDataScalar< ScalarT > terra::fv::hex::operators::FCTLimiterKernel< ScalarT >::antidiff_ |
\(\tilde{f}_{ij}\): pre-scaled antidiff fluxes.
|
staticconstexpr |
| grid::Grid4DDataScalar< ScalarT > terra::fv::hex::operators::FCTLimiterKernel< ScalarT >::R_minus_ |
Output: \(R_i^-\) correction factor.
| grid::Grid4DDataScalar< ScalarT > terra::fv::hex::operators::FCTLimiterKernel< ScalarT >::R_plus_ |
Output: \(R_i^+\) correction factor.
| grid::Grid4DDataScalar< ScalarT > terra::fv::hex::operators::FCTLimiterKernel< ScalarT >::T_L_ |
\(T^L\) with ghost layers filled.