4#include <AMReX_Config.H>
72 FAB const& operator[] (const
MFIter& mfi) const noexcept {
return m_mf[mfi]; }
88 [[nodiscard]]
Box fabbox (
int K)
const noexcept {
return m_mf.fabbox(K); }
90 [[nodiscard]]
int size () const noexcept {
return m_mf.size(); }
95 {
return m_mf.DistributionMap(); }
98 [[nodiscard]] MF
const&
multiFab () const noexcept {
return m_mf; }
100 [[nodiscard]]
int nComp () const noexcept {
return m_mf.nComp(); }
113 bool deterministic =
false);
115 void copyTo (MF& dest,
int ngrow,
int scomp,
int dcomp,
int ncomp,
118 void plusTo (MF& dest,
int ngrow,
int scomp,
int dcomp,
int ncomp,
127 int scomp,
int dcomp,
int ncomp);
132 int dcomp,
int ncomp,
int ngrow);
136 void write (
const std::string& name)
const;
139 void read (
const std::string& name);
152 template <
typename MF>
157template <
typename MF>
159 : m_mf(grids,dmap,ncomp,0)
162template <
typename MF>
166 m_mf.define(grids, dm, ncomp, 0);
169template <
typename MF>
175#pragma omp parallel if (Gpu::notInLaunchRegion())
178 const Box& bx = fsi.validbox();
179 auto const srcfab = src.
array(fsi);
180 auto dstfab = this->array(fsi);
183 dstfab(i,j,k,n+dcomp) = srcfab(i,j,k,n+scomp);
187 m_mf.ParallelCopy(src.
m_mf,scomp,dcomp,ncomp);
192template <
typename MF>
198 m_mf.ParallelCopy(src,scomp,dcomp,ncomp,ngrow,0,period);
202template <
typename MF>
208#pragma omp parallel if (Gpu::notInLaunchRegion())
211 const Box& bx = fsi.validbox();
212 auto const srcfab = src.
array(fsi);
213 auto dstfab = this->array(fsi);
216 dstfab(i,j,k,n+dcomp) += srcfab(i,j,k,n+scomp);
220 amrex::Abort(
"FabSetT<MF>::plusFrom: parallel plusFrom not supported");
225template <
typename MF>
231 m_mf.ParallelCopy(src,scomp,dcomp,ncomp,
IntVect(ngrow),
IntVect(0),period,
236template <
typename MF>
242 dest.ParallelCopy(m_mf,scomp,dcomp,ncomp,0,ngrow,period);
245template <
typename MF>
254template <
typename MF>
258 const int ncomp =
nComp();
260#pragma omp parallel if (Gpu::notInLaunchRegion())
263 const Box& bx = fsi.validbox();
264 auto fab = this->array(fsi);
272template <
typename MF>
277#pragma omp parallel if (Gpu::notInLaunchRegion())
280 const Box& bx = fsi.validbox();
281 auto fab = this->array(fsi);
284 fab(i,j,k,n+comp) = val;
291template <
typename MF>
294 int scomp,
int dcomp,
int ncomp)
299#pragma omp parallel if (Gpu::notInLaunchRegion())
303 const Box& bx = fsi.validbox();
304 auto const srcfab = src.
array(fsi);
305 auto dstfab = this->array(fsi);
308 dstfab(i,j,k,n+dcomp) = a*dstfab(i,j,k,n+dcomp) + b*srcfab(i,j,k,n+scomp);
316template <
typename MF>
320 int dcomp,
int ncomp,
int ngrow)
323 BL_ASSERT(mfa.nGrowVect().allGE(ngrow) && mfb.nGrowVect().allGE(ngrow));
324 BL_ASSERT(mfa.boxArray() == mfb.boxArray());
333#pragma omp parallel if (Gpu::notInLaunchRegion())
337 const Box& bx = mfi.validbox();
338 auto afab = bdrya.array(mfi);
339 auto bfab = bdryb.array(mfi);
342 afab(i,j,k,n) = huge;
343 bfab(i,j,k,n) = huge;
347 bdrya.ParallelCopy(mfa,a_comp,0,ncomp,ngrow,0);
348 bdryb.ParallelCopy(mfb,b_comp,0,ncomp,ngrow,0);
351#pragma omp parallel if (Gpu::notInLaunchRegion())
355 const Box& bx = fsi.validbox();
356 auto const afab = bdrya.array(fsi);
357 auto const bfab = bdryb.array(fsi);
358 auto dfab = this->array(fsi);
361 dfab(i,j,k,n+dcomp) = a*afab(i,j,k,n) + b*bfab(i,j,k,n);
368template <
typename MF>
379template <
typename MF>
389template <
typename MF>
395 int ncomp = dst.
nComp();
397#pragma omp parallel if (Gpu::notInLaunchRegion())
400 const Box& bx = fsi.validbox();
401 auto const srcfab = src.
array(fsi);
402 auto dstfab = dst.
array(fsi);
405 dstfab(i,j,k,n) = srcfab(i,j,k,n);
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:567
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
@ ADD
Definition AMReX_FabArrayBase.H:394
Definition AMReX_FabSet.H:150
FabSetIter(const FabSetT< MF > &fs)
Definition AMReX_FabSet.H:153
A FabSet is a group of FArrayBox's. The grouping is designed specifically to represent regions along ...
Definition AMReX_FabSet.H:46
Array4< value_type const > const_array(const MFIter &mfi) const noexcept
Definition AMReX_FabSet.H:81
MultiArray4< value_type const > arrays() const noexcept
Definition AMReX_FabSet.H:84
int nComp() const noexcept
Definition AMReX_FabSet.H:100
Array4< value_type const > const_array(int i) const noexcept
Definition AMReX_FabSet.H:82
MF & multiFab() noexcept
Definition AMReX_FabSet.H:97
Array4< value_type const > array(int i) const noexcept
Definition AMReX_FabSet.H:79
static void Copy(FabSetT< MF > &dst, const FabSetT< MF > &src)
Definition AMReX_FabSet.H:391
FabSetT< MF > & copyFrom(const FabSetT< MF > &src, int scomp, int dcomp, int ncomp)
Definition AMReX_FabSet.H:171
FabSetT< MF > & plusFrom(const FabSetT< MF > &src, int scomp, int dcomp, int ncomp)
Definition AMReX_FabSet.H:204
void read(const std::string &name)
Read (used for reading from checkpoint)
Definition AMReX_FabSet.H:381
void clear()
Definition AMReX_FabSet.H:102
FAB const & operator[](const MFIter &mfi) const noexcept
Definition AMReX_FabSet.H:72
void plusTo(MF &dest, int ngrow, int scomp, int dcomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic()) const
Definition AMReX_FabSet.H:247
MultiArray4< value_type > arrays() noexcept
Definition AMReX_FabSet.H:85
FabSetT< MF > & linComb(value_type a, value_type b, const FabSetT< MF > &src, int scomp, int dcomp, int ncomp)
Linear combination: this := a*mfa + b*mfb.
Definition AMReX_FabSet.H:293
Box fabbox(int K) const noexcept
Definition AMReX_FabSet.H:88
Array4< value_type > array(int i) noexcept
Definition AMReX_FabSet.H:80
MultiArray4< value_type const > const_arrays() const noexcept
Definition AMReX_FabSet.H:86
MF m_mf
Definition AMReX_FabSet.H:145
const BoxArray & boxArray() const noexcept
Definition AMReX_FabSet.H:92
void write(const std::string &name) const
Write (used for writing to checkpoint)
Definition AMReX_FabSet.H:370
FabSetT() noexcept=default
The default constructor – you must later call define().
typename FabDataType< MF >::value_type value_type
Definition AMReX_FabSet.H:50
void setVal(value_type val)
Definition AMReX_FabSet.H:256
int size() const noexcept
Definition AMReX_FabSet.H:90
Array4< value_type const > array(const MFIter &mfi) const noexcept
Definition AMReX_FabSet.H:77
typename FabDataType< MF >::fab_type FAB
Definition AMReX_FabSet.H:51
Array4< value_type > array(const MFIter &mfi) noexcept
Definition AMReX_FabSet.H:78
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_FabSet.H:94
void define(const BoxArray &grids, const DistributionMapping &dmap, int ncomp)
Define a FabSetT<MF> constructed via default constructor.
Definition AMReX_FabSet.H:164
MF const & multiFab() const noexcept
Definition AMReX_FabSet.H:98
void copyTo(MF &dest, int ngrow, int scomp, int dcomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic()) const
Definition AMReX_FabSet.H:238
Flux Register.
Definition AMReX_FluxRegister.H:20
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:63
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:147
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
static void Read(FabArray< FArrayBox > &mf, const std::string &name, const char *faHeader=nullptr, int coordinatorProc=ParallelDescriptor::IOProcessorNumber(), int allow_empty_mf=0)
Read a FabArray<FArrayBox> from disk written using VisMF::Write(). If the FabArray<FArrayBox> fafab h...
Definition AMReX_VisMF.cpp:1583
static Long Write(const FabArray< FArrayBox > &mf, const std::string &name, VisMF::How how=NFiles, bool set_ghost=false)
Write a FabArray<FArrayBox> to disk in a "smart" way. Returns the total number of bytes written on th...
Definition AMReX_VisMF.cpp:979
static void AsyncWrite(const FabArray< FArrayBox > &mf, const std::string &mf_name, bool valid_cells_only=false)
Definition AMReX_VisMF.cpp:2314
bool UseAsyncOut()
Definition AMReX_AsyncOut.cpp:70
Definition AMReX_Amr.cpp:49
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2854
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2869
bool match(const BoxArray &x, const BoxArray &y)
Note that two BoxArrays that match are not necessarily equal.
Definition AMReX_BoxArray.cpp:1927
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2864
Definition AMReX_Array4.H:61
Definition AMReX_FabDataType.H:9
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_FabArray.H:153