WarpX
Functions
CoarsenIO Namespace Reference

Functions

AMREX_GPU_DEVICE AMREX_FORCE_INLINE Real Interp (Array4< Real const > const &arr_src, GpuArray< int, 3 > const &sf, GpuArray< int, 3 > const &sc, GpuArray< int, 3 > const &cr, const int i, const int j, const int k, const int comp)
 Interpolates the floating point data contained in the source Array4 arr_src, extracted from a fine MultiFab, by averaging over either 1 point or 2 equally distant points. More...
 
void Loop (MultiFab &mf_dst, const MultiFab &mf_src, const int dcomp, const int scomp, const int ncomp, const int ngrow, const IntVect crse_ratio=IntVect(1))
 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 (MultiFab &mf_dst, const MultiFab &mf_src, const int dcomp, const int scomp, const int ncomp, const int ngrow, const IntVect crse_ratio=IntVect(1))
 Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in the fine MultiFab mf_src. More...
 

Function Documentation

◆ Coarsen()

void CoarsenIO::Coarsen ( MultiFab &  mf_dst,
const MultiFab &  mf_src,
const int  dcomp,
const int  scomp,
const int  ncomp,
const int  ngrow,
const IntVect  crse_ratio = IntVect(1) 
)

Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in the fine MultiFab mf_src.

Parameters
[in,out]mf_dstcoarsened MultiFab containing the floating point data to be filled by interpolating the fine MultiFab mf_src
[in]mf_srcfine MultiFab containing the floating point data to be interpolated
[in]dcompoffset for the fourth component of the coarsened Array4 object, extracted from its MultiFab, where the interpolated values will be stored
[in]scompoffset for the fourth component of the fine Array4 object, extracted from its MultiFab, containing the data to be interpolated
[in]ncompnumber of components to loop over for the coarsened Array4 extracted from the coarsened MultiFab mf_dst
[in]ngrownumber of guard cells to fill
[in]crse_ratiocoarsening ratio between the fine MultiFab mf_src and the coarsened MultiFab mf_dst along each spatial direction

◆ Interp()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE Real CoarsenIO::Interp ( Array4< Real const > const &  arr_src,
GpuArray< int, 3 > const &  sf,
GpuArray< int, 3 > const &  sc,
GpuArray< int, 3 > const &  cr,
const int  i,
const int  j,
const int  k,
const int  comp 
)

Interpolates the floating point data contained in the source Array4 arr_src, extracted from a fine MultiFab, by averaging over either 1 point or 2 equally distant points.

Parameters
[in]arr_srcfloating point data to be interpolated
[in]sfstaggering of the source fine MultiFab
[in]scstaggering of the destination coarsened MultiFab
[in]crcoarsening ratio along each spatial direction
[in]iindex along x of the coarsened Array4 to be filled
[in]jindex along y of the coarsened Array4 to be filled
[in]kindex along z of the coarsened Array4 to be filled
[in]compindex along the fourth component of the Array4 arr_src containing the data to be interpolated
Returns
interpolated field at cell (i,j,k) of a coarsened Array4

◆ Loop()

void CoarsenIO::Loop ( MultiFab &  mf_dst,
const MultiFab &  mf_src,
const int  dcomp,
const int  scomp,
const int  ncomp,
const int  ngrow,
const IntVect  crse_ratio = IntVect(1) 
)

Loops over the boxes of the coarsened MultiFab mf_dst and fills them by interpolating the data contained in the fine MultiFab mf_src.

Parameters
[in,out]mf_dstcoarsened MultiFab containing the floating point data to be filled by interpolating the source fine MultiFab
[in]mf_srcfine MultiFab containing the floating point data to be interpolated
[in]dcompoffset for the fourth component of the coarsened Array4 object, extracted from its MultiFab, where the interpolated values will be stored
[in]scompoffset for the fourth component of the fine Array4 object, extracted from its MultiFab, containing the data to be interpolated
[in]ncompnumber of components to loop over for the coarsened Array4 extracted from the coarsened MultiFab mf_dst
[in]ngrownumber of guard cells to fill
[in]crse_ratiocoarsening ratio between the fine MultiFab mf_src and the coarsened MultiFab mf_dst along each spatial direction