Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::EdgeFluxRegister Class Reference

#include <AMReX_EdgeFluxRegister.H>

Public Member Functions

 EdgeFluxRegister ()=default
 Construct an empty register; call define() before using it.
 
 EdgeFluxRegister (const BoxArray &fba, const BoxArray &cba, const DistributionMapping &fdm, const DistributionMapping &cdm, const Geometry &fgeom, const Geometry &cgeom, int nvar=1)
 Fully construct the register for a coarse/fine pair.
 
void define (const BoxArray &fba, const BoxArray &cba, const DistributionMapping &fdm, const DistributionMapping &cdm, const Geometry &fgeom, const Geometry &cgeom, int nvar=1)
 Define the register using fba, cba, fdm, cdm, fgeom, cgeom, and component count nvar.
 
void reset ()
 Reset accumulated state at the start of a coarse step.
 
void CrseAdd (MFIter const &mfi, const Array< FArrayBox const *, 3 > &E_crse, Real dt_crse)
 Accumulate coarse edge-centered fluxes for tile mfi.
 
void FineAdd (MFIter const &mfi, const Array< FArrayBox const *, 3 > &E_fine, Real dt_fine)
 Accumulate fine edge-centered fluxes for tile mfi.
 
void Reflux (Array< MultiFab *, 3 > const &B_crse) const
 Apply the accumulated edge fluxes to coarse face-centered data.
 

Detailed Description

Edge Flux Register for Constrained Transport

This Flux Register is useful for solving system like dB/dt + curl E = 0 on a staggered mesh. (Here d is of course partial derivation.) B is a vector on cell faces, and E is a vector on cell edges. In 2D, E has only one component, Ez, and it is on the nodes of a 2d mesh.

At the beginning of a coarse step, reset() is called. In MFIter for the coarse level advance, CrseAdd is called with coarse flux (i.e., E). The flux is not scaled. In MFIter for the fine level advance, FineAdd is called. After the fine level finishes its time steps, Reflux is called to update the coarse level B on the coarse/fine boundary. The user is also expected to call this version of average_down_faces from AMReX_MultiFabUtil.H to synchronize the coarse level data with the fine level.

template <typename MF, std::enable_if_t<IsFabArray<MF>::value,int>>
const IntVect& ratio, const Geometry& crse_geom)
Array4< Real > fine
Definition AMReX_InterpFaceRegister.cpp:90
Array4< Real const > crse
Definition AMReX_InterpFaceRegister.cpp:92
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
std::array< T, N > Array
Definition AMReX_Array.H:26
void average_down_faces(const Vector< const MF * > &fine, const Vector< MF * > &crse, const IntVect &ratio, int ngcrse=0)
Average fine face-based FabArray onto crse face-based FabArray.
Definition AMReX_MultiFabUtil.H:1033

Note that both CrseAdd and FineAdd are async for GPU builds. That means it's the user's responsibility to keep the FArrayBox arguments alive or call Gpu::streamSynchronize() when necessary.

Because staggered grids are used, tiling could be very confusing. To avoid confusion, this class assumes that tiling is not enabled for the MFIter loop containing calls to CrseAdd and FineAdd.

If the equation has an extra factor due to the choice of units, the factor can be absorbed into dt. If we have v x B instead of E, the sign can also been absorbed into dt. Note that whatever the choice of sign is, the dt arguments passed to CrseAdd and FineAdd should have the same sign.

We try to keep the interface simple by not providing overloads that specify the component index. If the user's data does not start with component 0, it can be worked around by creating alias FArrayBox and MultiFab.

Constructor & Destructor Documentation

◆ EdgeFluxRegister() [1/2]

amrex::EdgeFluxRegister::EdgeFluxRegister ( )
default

Construct an empty register; call define() before using it.

◆ EdgeFluxRegister() [2/2]

amrex::EdgeFluxRegister::EdgeFluxRegister ( const BoxArray fba,
const BoxArray cba,
const DistributionMapping fdm,
const DistributionMapping cdm,
const Geometry fgeom,
const Geometry cgeom,
int  nvar = 1 
)

Fully construct the register for a coarse/fine pair.

Parameters
fbaFine-level BoxArray
cbaCoarse-level BoxArray
fdmDistribution mapping for the fine grids.
cdmDistribution mapping for the coarse grids.
fgeomFine-level Geometry.
cgeomCoarse-level Geometry.
nvarNumber of edge-centered components.

Member Function Documentation

◆ CrseAdd()

void amrex::EdgeFluxRegister::CrseAdd ( MFIter const &  mfi,
const Array< FArrayBox const *, 3 > &  E_crse,
Real  dt_crse 
)

Accumulate coarse edge-centered fluxes for tile mfi.

Parameters
mfiIterator identifying the coarse tile.
E_crseArray of pointers (Ex,Ey,Ez) to coarse-time edge data.
dt_crseTime step scaling associated with the coarse advance.

◆ define()

void amrex::EdgeFluxRegister::define ( const BoxArray fba,
const BoxArray cba,
const DistributionMapping fdm,
const DistributionMapping cdm,
const Geometry fgeom,
const Geometry cgeom,
int  nvar = 1 
)

Define the register using fba, cba, fdm, cdm, fgeom, cgeom, and component count nvar.

◆ FineAdd()

void amrex::EdgeFluxRegister::FineAdd ( MFIter const &  mfi,
const Array< FArrayBox const *, 3 > &  E_fine,
Real  dt_fine 
)

Accumulate fine edge-centered fluxes for tile mfi.

Parameters
mfiIterator identifying the fine tile.
E_fineArray of pointers (Ex,Ey,Ez) to fine-time edge data.
dt_fineTime step scaling associated with the fine advance.

◆ Reflux()

void amrex::EdgeFluxRegister::Reflux ( Array< MultiFab *, 3 > const &  B_crse) const

Apply the accumulated edge fluxes to coarse face-centered data.

Parameters
B_crseArray of coarse face MultiFabs to update (Bx, By, Bz).

◆ reset()

void amrex::EdgeFluxRegister::reset ( )

Reset accumulated state at the start of a coarse step.


The documentation for this class was generated from the following files: