1#ifndef AMREX_MULTICUTFAB_H_
2#define AMREX_MULTICUTFAB_H_
3#include <AMReX_Config.H>
29 bool shared=
false,
Arena* ar =
nullptr)
33 :
FArrayBox(rhs, make_type, scomp, ncomp) {}
45 return copyFromMem<run_on, BUF>(
box(), 0,
nComp(), src);
48 template <RunOn run_on,
typename BUF = CutFab::value_type>
54 if (
dptr !=
nullptr) {
55 return FArrayBox::copyFromMem<run_on, BUF>(dstbox, dstcomp, numcomp, src);
57 return sizeof(BUF)*
static_cast<std::size_t
>(dstbox.
numPts()*numcomp);
61 template <RunOn run_on>
69 if (
dptr !=
nullptr) {
70 FArrayBox::copy<run_on>(src,srcbox,srccomp,destbox,destcomp,numcomp);
75 template <RunOn run_on>
78 if (
dptr !=
nullptr) {
79 FArrayBox::copy<run_on>(src, bx, scomp, dcomp, ncomp);
120 return m_data.arrays();
124 return m_data.const_arrays();
128 return m_data.const_arrays();
132 bool ok (
const MFIter& mfi)
const noexcept;
135 bool ok (
int global_box_index)
const noexcept;
145 int nComp () const noexcept {
return m_data.nComp(); }
146 int nGrow () const noexcept {
return m_data.nGrow(); }
149 void ParallelCopy (
const MultiCutFab& src,
int scomp,
int dcomp,
int ncomp,
int sng,
int dng,
A virtual base class for objects that manage their own dynamic memory allocation.
Definition AMReX_Arena.H:132
const Box & box() const noexcept
Returns the domain (box) where the array is defined.
Definition AMReX_BaseFab.H:294
Real value_type
Definition AMReX_BaseFab.H:195
Real * dptr
The data pointer.
Definition AMReX_BaseFab.H:1164
int nComp() const noexcept
Returns the number of components.
Definition AMReX_BaseFab.H:280
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
__host__ __device__ Long numPts() const noexcept
Return the number of points contained in the BoxND.
Definition AMReX_Box.H:356
FArrayBox specialization used to store cut-cell data (may skip allocation).
Definition AMReX_MultiCutFab.H:18
CutFab & copy(const CutFab &src, const Box &bx, SrcComp scomp, DestComp dcomp, NumComps ncomp)
Definition AMReX_MultiCutFab.H:76
CutFab(Arena *ar)
Definition AMReX_MultiCutFab.H:23
~CutFab() noexcept override=default
CutFab(const Box &b, int n, Arena *ar)
Definition AMReX_MultiCutFab.H:25
CutFab(CutFab const &rhs, MakeType make_type, int scomp, int ncomp)
Definition AMReX_MultiCutFab.H:32
CutFab(const Box &b, int ncomps, bool alloc=true, bool shared=false, Arena *ar=nullptr)
Definition AMReX_MultiCutFab.H:28
std::size_t copyFromMem(const Box &dstbox, int dstcomp, int numcomp, const void *src)
Definition AMReX_MultiCutFab.H:49
CutFab & copy(const CutFab &src, const Box &srcbox, int srccomp, const Box &destbox, int destcomp, int numcomp)
Definition AMReX_MultiCutFab.H:62
std::size_t copyFromMem(const void *src)
Definition AMReX_MultiCutFab.H:44
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:350
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
FabArray façade that only allocates cut-cell Fabs (skipping regular cells).
Definition AMReX_MultiCutFab.H:87
MultiCutFab(const MultiCutFab &rhs)=delete
MultiArray4< Real > arrays() noexcept
Definition AMReX_MultiCutFab.H:119
MultiArray4< Real const > arrays() const noexcept
Definition AMReX_MultiCutFab.H:123
void setVal(Real val)
Fill every stored CutFab with a constant val.
Definition AMReX_MultiCutFab.cpp:131
Array4< Real > array(const MFIter &mfi) noexcept
Array4 accessors (mutable and const) matching MFIter order.
Definition AMReX_MultiCutFab.cpp:112
const CutFab & operator[](const MFIter &mfi) const noexcept
Access the CutFab for mfi.
Definition AMReX_MultiCutFab.cpp:70
int nComp() const noexcept
Definition AMReX_MultiCutFab.H:145
MultiCutFab(MultiCutFab &&rhs) noexcept=default
const FabArray< CutFab > & data() const noexcept
Definition AMReX_MultiCutFab.H:141
FabArray< CutFab > & data() noexcept
Definition AMReX_MultiCutFab.H:140
void define(const BoxArray &ba, const DistributionMapping &dm, int ncomp, int ngrow, const FabArray< EBCellFlagFab > &cellflags)
Definition AMReX_MultiCutFab.cpp:62
MultiArray4< Real const > const_arrays() const noexcept
Definition AMReX_MultiCutFab.H:127
MultiFab ToMultiFab(Real regular_value, Real covered_value) const
Materialize the cut-cell data into a dense MultiFab.
Definition AMReX_MultiCutFab.cpp:156
Array4< Real const > const_array(const MFIter &mfi) const noexcept
Definition AMReX_MultiCutFab.cpp:98
const BoxArray & boxArray() const noexcept
Definition AMReX_MultiCutFab.H:143
bool ok(const MFIter &mfi) const noexcept
Is it OK to call operator[] with this MFIter?
Definition AMReX_MultiCutFab.cpp:119
int nGrow() const noexcept
Definition AMReX_MultiCutFab.H:146
MultiCutFab & operator=(const MultiCutFab &rhs)=delete
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_MultiCutFab.H:144
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
static const Periodicity & NonPeriodic() noexcept
Definition AMReX_Periodicity.cpp:52
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
Definition AMReX_Amr.cpp:49
MakeType
Definition AMReX_MakeType.H:7
RunOn
Definition AMReX_GpuControl.H:69
A multidimensional array accessor.
Definition AMReX_Array4.H:283
void * alloc(std::size_t sz) const noexcept
Definition AMReX_DataAllocator.H:16
Definition AMReX_BaseFab.H:77
Definition AMReX_FabArray.H:154
Definition AMReX_BaseFab.H:83
Definition AMReX_BaseFab.H:71