WarpX
|
Functions | |
void | Loop (amrex::MultiFab &mf_dst, amrex::MultiFab const &mf_src, int ncomp, amrex::IntVect ngrow, amrex::IntVect crse_ratio) |
Loops over the boxes of the coarsened MultiFab mf_dst and fills them by interpolating the data contained in the fine MultiFab mf_src . More... | |
void | Coarsen (amrex::MultiFab &mf_dst, amrex::MultiFab const &mf_src, amrex::IntVect crse_ratio) |
Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in the fine MultiFab mf_src . More... | |
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real | Interp (amrex::Array4< amrex::Real const > const &arr_src, amrex::GpuArray< int, 3 > const &sf, amrex::GpuArray< int, 3 > const &sc, amrex::GpuArray< int, 3 > const &cr, int const i, int const j, int const k, int const comp) |
Interpolates the floating point data contained in the source Array4 arr_src , extracted from a fine MultiFab, with weights defined in such a way that the total charge is preserved. More... | |
Mesh Coarsening by Averaging
These methods are mostly used for mesh-refinement.
void ablastr::coarsen::average::Coarsen | ( | amrex::MultiFab & | mf_dst, |
amrex::MultiFab const & | mf_src, | ||
amrex::IntVect | crse_ratio | ||
) |
Stores in the coarsened MultiFab mf_dst
the values obtained by interpolating the data contained in the fine MultiFab mf_src
.
[in,out] | mf_dst | coarsened MultiFab containing the floating point data to be filled by interpolating the fine MultiFab mf_src |
[in] | mf_src | fine MultiFab containing the floating point data to be interpolated |
[in] | crse_ratio | coarsening ratio between the fine MultiFab mf_src and the coarsened MultiFab mf_dst along each spatial direction |
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ablastr::coarsen::average::Interp | ( | amrex::Array4< amrex::Real const > const & | arr_src, |
amrex::GpuArray< int, 3 > const & | sf, | ||
amrex::GpuArray< int, 3 > const & | sc, | ||
amrex::GpuArray< int, 3 > const & | cr, | ||
int const | i, | ||
int const | j, | ||
int const | k, | ||
int const | comp | ||
) |
Interpolates the floating point data contained in the source Array4 arr_src
, extracted from a fine MultiFab, with weights defined in such a way that the total charge is preserved.
The input (sf) and output (sc) staggering need to be the same.
[in] | arr_src | floating point data to be interpolated |
[in] | sf | staggering of the source fine MultiFab |
[in] | sc | staggering of the destination coarsened MultiFab |
[in] | cr | coarsening ratio along each spatial direction |
[in] | i | index along x of the coarsened Array4 to be filled |
[in] | j | index along y of the coarsened Array4 to be filled |
[in] | k | index along z of the coarsened Array4 to be filled |
[in] | comp | index along the fourth component of the Array4 arr_src containing the data to be interpolated |
void ablastr::coarsen::average::Loop | ( | amrex::MultiFab & | mf_dst, |
amrex::MultiFab const & | mf_src, | ||
int | ncomp, | ||
amrex::IntVect | ngrow, | ||
amrex::IntVect | crse_ratio | ||
) |
Loops over the boxes of the coarsened MultiFab mf_dst
and fills them by interpolating the data contained in the fine MultiFab mf_src
.
[in,out] | mf_dst | coarsened MultiFab containing the floating point data to be filled by interpolating the source fine MultiFab |
[in] | mf_src | fine MultiFab containing the floating point data to be interpolated |
[in] | ncomp | number of components to loop over for the coarsened Array4 extracted from the coarsened MultiFab mf_dst |
[in] | ngrow | number of guard cells to fill along each spatial direction |
[in] | crse_ratio | coarsening ratio between the fine MultiFab mf_src and the coarsened MultiFab mf_dst along each spatial direction |