![]() |
Block-Structured AMR Software Framework
|
#include <AMReX_YAFluxRegister.H>
Public Types | |
| enum | CellType : int { crse_cell = 0 , crse_fine_boundary_cell , fine_cell } |
| using | T = typename MF::value_type |
| using | FAB = typename MF::fab_type |
Public Member Functions | |
| YAFluxRegisterT ()=default | |
| Construct an empty register; call define() before use. | |
| YAFluxRegisterT (const BoxArray &fba, const BoxArray &cba, const DistributionMapping &fdm, const DistributionMapping &cdm, const Geometry &fgeom, const Geometry &cgeom, const IntVect &ref_ratio, int fine_lev, int nvar) | |
| Fully construct the register for a fine level. | |
| void | define (const BoxArray &fba, const BoxArray &cba, const DistributionMapping &fdm, const DistributionMapping &cdm, const Geometry &fgeom, const Geometry &cgeom, const IntVect &ref_ratio, int fine_lev, int nvar) |
| Define the register using the same arguments as the constructor. | |
| void | reset () |
| Reset stored fluxes and flags before a coarse advance. | |
| void | CrseAdd (const MFIter &mfi, const std::array< FAB const *, 3 > &flux, const Real *dx, Real dt, RunOn runon) noexcept |
Add coarse-level fluxes for the tile indicated by mfi. | |
| void | CrseAdd (const MFIter &mfi, const std::array< FAB const *, 3 > &flux, const Real *dx, Real dt, int srccomp, int destcomp, int numcomp, RunOn runon) noexcept |
| Add coarse-level fluxes with component slicing. | |
| void | FineAdd (const MFIter &mfi, const std::array< FAB const *, 3 > &flux, const Real *dx, Real dt, RunOn runon) noexcept |
Add fine-level fluxes for the tile identified by mfi. | |
| void | FineAdd (const MFIter &mfi, const std::array< FAB const *, 3 > &a_flux, const Real *dx, Real dt, int srccomp, int destcomp, int numcomp, RunOn runon) noexcept |
| Add fine-level fluxes with component slicing. | |
| void | Reflux (MF &state, int dc=0) |
Apply the accumulated flux divergence to the coarse MultiFab state. | |
| void | Reflux (MF &state, int srccomp, int destcomp, int numcomp) |
| Component-aware version of Reflux writing components [destcomp,destcomp+numcomp). | |
| bool | CrseHasWork (const MFIter &mfi) const noexcept |
Return true if coarse tile mfi abuts a coarse/fine boundary. | |
| bool | FineHasWork (const MFIter &mfi) const noexcept |
Return true if fine tile mfi has flux registers to update. | |
| MF & | getFineData () |
| Mutable access to fine-side accumulation data. | |
| MF & | getCrseData () |
| Mutable access to coarse-side accumulation data. | |
| void | setCrseVolume (MF const *cvol) |
| void | setDeterministic (bool flag) |
Enable deterministic mode for GPU operations via flag (uses deterministic reductions). | |
| bool | getDeterministic () const |
Protected Attributes | |
| MF | m_crse_data |
| iMultiFab | m_crse_flag |
| Vector< int > | m_crse_fab_flag |
| MF | m_cfpatch |
| This is built on crse/fine patches. | |
| MF | m_cfp_mask |
| Vector< Vector< FAB * > > | m_cfp_fab |
| The size of this is (# of local fine grids (# of crse/fine patches for that grid)) | |
| Vector< int > | m_cfp_localindex |
| Geometry | m_fine_geom |
| Geometry | m_crse_geom |
| IntVect | m_ratio |
| int | m_fine_level |
| int | m_ncomp |
| MF const * | m_cvol = nullptr |
| bool | m_deterministic = false |
YAFluxRegister is yet another Flux Register for refluxing.
At the beginning of a coarse step, reset() is called. In MFIter for the coarse level advance, CrseAdd is called with coarse flux. The flux is not scaled. In MFIter for the fine level advance, FineAdd is called. After the fine level finished its time steps, Reflux is called to update the coarse cells next to the coarse/fine boundary.
| using amrex::YAFluxRegisterT< MF >::FAB = typename MF::fab_type |
| using amrex::YAFluxRegisterT< MF >::T = typename MF::value_type |
| enum amrex::YAFluxRegisterT::CellType : int |
|
default |
Construct an empty register; call define() before use.
| amrex::YAFluxRegisterT< MF >::YAFluxRegisterT | ( | const BoxArray & | fba, |
| const BoxArray & | cba, | ||
| const DistributionMapping & | fdm, | ||
| const DistributionMapping & | cdm, | ||
| const Geometry & | fgeom, | ||
| const Geometry & | cgeom, | ||
| const IntVect & | ref_ratio, | ||
| int | fine_lev, | ||
| int | nvar | ||
| ) |
Fully construct the register for a fine level.
|
noexcept |
Add coarse-level fluxes with component slicing.
| mfi | Iterator describing the coarse tile. |
| flux | Array of per-direction flux FABs. |
| dx | Cell spacing array. |
| dt | Coarse timestep. |
| srccomp | Starting component in the fluxes. |
| destcomp | Destination component in the register. |
| numcomp | Number of components. |
| runon | Execution context (host/device). |
|
noexcept |
Add coarse-level fluxes for the tile indicated by mfi.
| mfi | Iterator describing the coarse tile. |
| flux | Array of per-direction flux FABs. |
| dx | Cell spacing array. |
| dt | Coarse timestep (with sign). |
| runon | Execution context (host/device). |
|
inlinenoexcept |
Return true if coarse tile mfi abuts a coarse/fine boundary.
| void amrex::YAFluxRegisterT< MF >::define | ( | const BoxArray & | fba, |
| const BoxArray & | cba, | ||
| const DistributionMapping & | fdm, | ||
| const DistributionMapping & | cdm, | ||
| const Geometry & | fgeom, | ||
| const Geometry & | cgeom, | ||
| const IntVect & | ref_ratio, | ||
| int | fine_lev, | ||
| int | nvar | ||
| ) |
Define the register using the same arguments as the constructor.
|
noexcept |
Add fine-level fluxes with component slicing.
| mfi | Iterator describing the fine tile. |
| a_flux | Array of per-direction flux FABs. |
| dx | Fine cell spacing array. |
| dt | Fine timestep. |
| srccomp | Starting component in the fluxes. |
| destcomp | Destination component in the register. |
| numcomp | Number of components. |
| runon | Execution context (host/device). |
|
noexcept |
Add fine-level fluxes for the tile identified by mfi.
| mfi | Iterator describing the fine tile. |
| flux | Array of per-direction flux FABs. |
| dx | Fine cell spacing array. |
| dt | Fine timestep (per substep). |
| runon | Execution context (host/device). |
|
inlinenoexcept |
Return true if fine tile mfi has flux registers to update.
| MF & amrex::YAFluxRegisterT< MF >::getCrseData | ( | ) |
Mutable access to coarse-side accumulation data.
|
inline |
| MF & amrex::YAFluxRegisterT< MF >::getFineData | ( | ) |
Mutable access to fine-side accumulation data.
| void amrex::YAFluxRegisterT< MF >::Reflux | ( | MF & | state, |
| int | dc = 0 |
||
| ) |
| void amrex::YAFluxRegisterT< MF >::Reflux | ( | MF & | state, |
| int | srccomp, | ||
| int | destcomp, | ||
| int | numcomp | ||
| ) |
Component-aware version of Reflux writing components [destcomp,destcomp+numcomp).
| void amrex::YAFluxRegisterT< MF >::reset | ( | ) |
Reset stored fluxes and flags before a coarse advance.
|
inline |
For curvilinear coordinates only, where the supplied flux has already been multiplied by area. Provide the coarse-volume MultiFab cvol and keep it alive for the register lifetime (YAFluxRegister does not copy it).
|
inline |
Enable deterministic mode for GPU operations via flag (uses deterministic reductions).
|
protected |
The size of this is (# of local fine grids (# of crse/fine patches for that grid))
|
protected |
|
protected |
|
protected |
This is built on crse/fine patches.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |