Kokkos kernel for the explicit low-order predictor and antidiffusive flux assembly. More...
#include <fct_advection_diffusion.hpp>
Public Types | |
| using | ScalarType = ScalarT |
| using | Vec3 = dense::Vec< ScalarT, 3 > |
| 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::Grid3DDataVec< ScalarT, 3 > | grid_ |
| grid::Grid2DDataScalar< ScalarT > | radii_ |
| grid::Grid4DDataVec< ScalarT, 3 > | cell_centers_ |
| grid::Grid4DDataVec< ScalarT, 3 > | vel_grid_ |
| grid::Grid4DDataScalar< ScalarT > | T_old_ |
| \(T^n\): scalar at time level \(n\) (ghost layers must be filled). | |
| grid::Grid4DDataScalar< ScalarT > | T_L_ |
| \(T^L\): low-order predictor output. | |
| grid::Grid5DDataScalar< ScalarT > | antidiff_ |
| \(\tilde{f}_{ij}\): pre-scaled antidiff flux per face, shape \([\ldots, 6]\). | |
| ScalarT | dt_ |
| Time step \(\Delta t\). | |
| ScalarT | diffusivity_ |
| Physical diffusivity \(\kappa \geq 0\); set to 0 for pure advection. | |
| bool | subtract_divergence_ = true |
| Whether to subtract the discrete divergence error \(T_i\,\nabla\cdot\mathbf{u}\). | |
Optional volumetric source term \f$f\f$ [T/time] | |
Source values per cell. A null (default-constructed) view means no source. The explicit predictor adds \(\Delta t \cdot f_i\) to \(T_i^L\) after the upwind/diffusion update. Physical units: same as \(T\) per unit time. | |
| grid::Grid4DDataScalar< ScalarT > | source_ |
| bool | has_source_ = false |
Static Public Attributes | |
| static constexpr int | num_neighbors = GH::num_neighbors |
Kokkos kernel for the explicit low-order predictor and antidiffusive flux assembly.
In a single pass over each cell \(i\) this kernel:
GeometryHelper::compute_geometry to obtain the face-normal velocity fluxes \(\beta_{ij}\), the cell volume \(M_{ii}\), and the area-weighted face normals \(\mathbf{S}_f^{(j)}\) for all 6 neighbours.\[ T_i^L = T_i^n - \frac{\Delta t}{M_{ii}} \left[ \beta_{ii}^+\,T_i^n + \sum_j \beta_{ij}^-\,T_j^n + \sum_j \kappa\, \frac{|\mathbf{S}_f^{(j)}|^2} {(\mathbf{x}_j-\mathbf{x}_i)\cdot\mathbf{S}_f^{(j)}} (T_i^n - T_j^n) \right], \]
where \(\beta^+ = \max(\beta,0)\) and \(\beta^- = \min(\beta,0)\).\[ \tilde{f}_{ij} = \frac{\Delta t}{M_{ii}}\,\frac{|\beta_{ij}|}{2}\,(T_i^n - T_j^n). \]
Physical diffusion is not included in \(\tilde{f}_{ij}\) — only purely advective antidiffusion enters the Zalesak limiter.cell_centers_) must be filled before launch. | ScalarT | Floating-point scalar type. |
| using terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::GH = fct_detail::GeometryHelper< ScalarT > |
| using terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::ScalarType = ScalarT |
| using terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::Vec3 = dense::Vec< ScalarT, 3 > |
|
inline |
| grid::Grid5DDataScalar< ScalarT > terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::antidiff_ |
\(\tilde{f}_{ij}\): pre-scaled antidiff flux per face, shape \([\ldots, 6]\).
| grid::Grid4DDataVec< ScalarT, 3 > terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::cell_centers_ |
| ScalarT terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::diffusivity_ |
Physical diffusivity \(\kappa \geq 0\); set to 0 for pure advection.
| ScalarT terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::dt_ |
Time step \(\Delta t\).
| grid::Grid3DDataVec< ScalarT, 3 > terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::grid_ |
| bool terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::has_source_ = false |
|
staticconstexpr |
| grid::Grid2DDataScalar< ScalarT > terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::radii_ |
| grid::Grid4DDataScalar< ScalarT > terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::source_ |
| bool terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::subtract_divergence_ = true |
Whether to subtract the discrete divergence error \(T_i\,\nabla\cdot\mathbf{u}\).
When true (default), the correction \(-T_i\,(\sum_j \beta_{ij})/M_{ii}\) is applied, converting the conservative form \(\nabla\cdot(\mathbf{u}T)\) into the advective form \(\mathbf{u}\cdot\nabla T\). This is the physically correct equation for temperature. When \(\nabla\cdot\mathbf{u} = 0\) the correction is exactly zero and costs nothing.
| grid::Grid4DDataScalar< ScalarT > terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::T_L_ |
\(T^L\): low-order predictor output.
| grid::Grid4DDataScalar< ScalarT > terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::T_old_ |
\(T^n\): scalar at time level \(n\) (ghost layers must be filled).
| grid::Grid4DDataVec< ScalarT, 3 > terra::fv::hex::operators::FCTPredictorKernel< ScalarT >::vel_grid_ |