Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_FabSet.H
Go to the documentation of this file.
1
2#ifndef AMREX_FABSET_H_
3#define AMREX_FABSET_H_
4#include <AMReX_Config.H>
5
6#include <AMReX_FabDataType.H>
7#include <AMReX_MultiFab.H>
9#include <AMReX_BLProfiler.H>
10#include <AMReX_VisMF.H>
11
12#ifdef AMREX_USE_OMP
13#include <omp.h>
14#endif
15
16#include <limits>
17
18namespace amrex {
19
44template <typename MF>
46{
47 friend class FabSetIter;
48 friend class FluxRegister;
49public:
52
53 //
55 FabSetT () noexcept = default;
56 //
58 FabSetT (const BoxArray& grids, const DistributionMapping& dmap, int ncomp);
59
60 ~FabSetT () = default;
61
62 FabSetT (FabSetT<MF>&& rhs) noexcept = default;
63
64 FabSetT (const FabSetT<MF>& rhs) = delete;
65 FabSetT<MF>& operator= (const FabSetT<MF>& rhs) = delete;
66 FabSetT<MF>& operator= (FabSetT<MF>&& rhs) = delete;
67
68 //
70 void define (const BoxArray& grids, const DistributionMapping& dmap, int ncomp);
71
72 FAB const& operator[] (const MFIter& mfi) const noexcept { return m_mf[mfi]; }
73 FAB & operator[] (const MFIter& mfi) noexcept { return m_mf[mfi]; }
74 FAB const& operator[] (int i) const noexcept { return m_mf[i]; }
75 FAB & operator[] (int i) noexcept { return m_mf[i]; }
76
77 Array4<value_type const> array (const MFIter& mfi) const noexcept { return m_mf.const_array(mfi); }
78 Array4<value_type > array (const MFIter& mfi) noexcept { return m_mf.array(mfi); }
79 Array4<value_type const> array (int i) const noexcept { return m_mf.const_array(i); }
80 Array4<value_type > array (int i) noexcept { return m_mf.array(i); }
81 Array4<value_type const> const_array (const MFIter& mfi) const noexcept { return m_mf.const_array(mfi); }
82 Array4<value_type const> const_array (int i) const noexcept { return m_mf.const_array(i); }
83
84 MultiArray4<value_type const> arrays () const noexcept { return m_mf.const_arrays(); }
85 MultiArray4<value_type > arrays () noexcept { return m_mf.arrays(); }
86 MultiArray4<value_type const> const_arrays () const noexcept { return m_mf.const_arrays(); }
87
88 [[nodiscard]] Box fabbox (int K) const noexcept { return m_mf.fabbox(K); }
89
90 [[nodiscard]] int size () const noexcept { return m_mf.size(); }
91
92 [[nodiscard]] const BoxArray& boxArray () const noexcept { return m_mf.boxArray(); }
93
94 [[nodiscard]] const DistributionMapping& DistributionMap () const noexcept
95 { return m_mf.DistributionMap(); }
96
97 [[nodiscard]] MF & multiFab () noexcept { return m_mf; }
98 [[nodiscard]] MF const& multiFab () const noexcept { return m_mf; }
99
100 [[nodiscard]] int nComp () const noexcept { return m_mf.nComp(); }
101
102 void clear () { m_mf.clear(); }
103
104 FabSetT<MF>& copyFrom (const FabSetT<MF>& src, int scomp, int dcomp, int ncomp);
105
106 FabSetT<MF>& copyFrom (const MF& src, int ngrow, int scomp, int dcomp, int ncomp,
107 const Periodicity& period = Periodicity::NonPeriodic());
108
109 FabSetT<MF>& plusFrom (const FabSetT<MF>& src, int scomp, int dcomp, int ncomp);
110
111 FabSetT<MF>& plusFrom (const MF& src, int ngrow, int scomp, int dcomp, int ncomp,
112 const Periodicity& period = Periodicity::NonPeriodic(),
113 bool deterministic = false);
114
115 void copyTo (MF& dest, int ngrow, int scomp, int dcomp, int ncomp,
116 const Periodicity& period = Periodicity::NonPeriodic()) const;
117
118 void plusTo (MF& dest, int ngrow, int scomp, int dcomp, int ncomp,
119 const Periodicity& period = Periodicity::NonPeriodic()) const;
120
121 void setVal (value_type val);
122
123 void setVal (value_type val, int comp, int num_comp);
124
127 int scomp, int dcomp, int ncomp);
128
130 FabSetT<MF>& linComb (value_type a, const MF& mfa, int a_comp,
131 value_type b, const MF& mfb, int b_comp,
132 int dcomp, int ncomp, int ngrow);
133
134 //
136 void write (const std::string& name) const;
137 //
139 void read (const std::string& name);
140
142 static void Copy (FabSetT<MF>& dst, const FabSetT<MF>& src);
143
144private:
146};
147
149 : public MFIter
150{
151public:
152 template <typename MF>
153 explicit FabSetIter (const FabSetT<MF>& fs)
154 : MFIter(fs.m_mf) { }
155};
156
157template <typename MF>
158FabSetT<MF>::FabSetT (const BoxArray& grids, const DistributionMapping& dmap, int ncomp)
159 : m_mf(grids,dmap,ncomp,0)
160{}
161
162template <typename MF>
163void
164FabSetT<MF>::define (const BoxArray& grids, const DistributionMapping& dm, int ncomp)
165{
166 m_mf.define(grids, dm, ncomp, 0);
167}
168
169template <typename MF>
171FabSetT<MF>::copyFrom (const FabSetT<MF>& src, int scomp, int dcomp, int ncomp)
172{
173 if (boxArray() == src.boxArray() && DistributionMap() == src.DistributionMap()) {
174#ifdef AMREX_USE_OMP
175#pragma omp parallel if (Gpu::notInLaunchRegion())
176#endif
177 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi) {
178 const Box& bx = fsi.validbox();
179 auto const srcfab = src.array(fsi);
180 auto dstfab = this->array(fsi);
181 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
182 {
183 dstfab(i,j,k,n+dcomp) = srcfab(i,j,k,n+scomp);
184 });
185 }
186 } else {
187 m_mf.ParallelCopy(src.m_mf,scomp,dcomp,ncomp);
188 }
189 return *this;
190}
191
192template <typename MF>
194FabSetT<MF>::copyFrom (const MF& src, int ngrow, int scomp, int dcomp, int ncomp,
195 const Periodicity& period)
196{
197 BL_ASSERT(boxArray() != src.boxArray());
198 m_mf.ParallelCopy(src,scomp,dcomp,ncomp,ngrow,0,period);
199 return *this;
200}
201
202template <typename MF>
204FabSetT<MF>::plusFrom (const FabSetT<MF>& src, int scomp, int dcomp, int ncomp)
205{
206 if (boxArray() == src.boxArray() && DistributionMap() == src.DistributionMap()) {
207#ifdef AMREX_USE_OMP
208#pragma omp parallel if (Gpu::notInLaunchRegion())
209#endif
210 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi) {
211 const Box& bx = fsi.validbox();
212 auto const srcfab = src.array(fsi);
213 auto dstfab = this->array(fsi);
214 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
215 {
216 dstfab(i,j,k,n+dcomp) += srcfab(i,j,k,n+scomp);
217 });
218 }
219 } else {
220 amrex::Abort("FabSetT<MF>::plusFrom: parallel plusFrom not supported");
221 }
222 return *this;
223}
224
225template <typename MF>
227FabSetT<MF>::plusFrom (const MF& src, int ngrow, int scomp, int dcomp, int ncomp,
228 const Periodicity& period, bool deterministic)
229{
230 BL_ASSERT(boxArray() != src.boxArray());
231 m_mf.ParallelCopy(src,scomp,dcomp,ncomp,IntVect(ngrow),IntVect(0),period,
232 FabArrayBase::ADD,nullptr,deterministic);
233 return *this;
234}
235
236template <typename MF>
237void
238FabSetT<MF>::copyTo (MF& dest, int ngrow, int scomp, int dcomp, int ncomp,
239 const Periodicity& period) const
240{
241 BL_ASSERT(boxArray() != dest.boxArray());
242 dest.ParallelCopy(m_mf,scomp,dcomp,ncomp,0,ngrow,period);
243}
244
245template <typename MF>
246void
247FabSetT<MF>::plusTo (MF& dest, int ngrow, int scomp, int dcomp, int ncomp,
248 const Periodicity& period) const
249{
250 BL_ASSERT(boxArray() != dest.boxArray());
251 dest.ParallelCopy(m_mf,scomp,dcomp,ncomp,0,ngrow,period,FabArrayBase::ADD);
252}
253
254template <typename MF>
255void
257{
258 const int ncomp = nComp();
259#ifdef AMREX_USE_OMP
260#pragma omp parallel if (Gpu::notInLaunchRegion())
261#endif
262 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi) {
263 const Box& bx = fsi.validbox();
264 auto fab = this->array(fsi);
265 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
266 {
267 fab(i,j,k,n) = val;
268 });
269 }
270}
271
272template <typename MF>
273void
274FabSetT<MF>::setVal (value_type val, int comp, int num_comp)
275{
276#ifdef AMREX_USE_OMP
277#pragma omp parallel if (Gpu::notInLaunchRegion())
278#endif
279 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi) {
280 const Box& bx = fsi.validbox();
281 auto fab = this->array(fsi);
282 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, num_comp, i, j, k, n,
283 {
284 fab(i,j,k,n+comp) = val;
285 });
286 }
287}
288
289// Linear combination this := a*this + b*src
290// Note: corresponding fabsets must be commensurate.
291template <typename MF>
294 int scomp, int dcomp, int ncomp)
295{
296 BL_ASSERT(size() == src.size());
297
298#ifdef AMREX_USE_OMP
299#pragma omp parallel if (Gpu::notInLaunchRegion())
300#endif
301 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi)
302 {
303 const Box& bx = fsi.validbox();
304 auto const srcfab = src.array(fsi);
305 auto dstfab = this->array(fsi);
306 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
307 {
308 dstfab(i,j,k,n+dcomp) = a*dstfab(i,j,k,n+dcomp) + b*srcfab(i,j,k,n+scomp);
309 });
310 }
311 return *this;
312}
313
314// Linear combination: this := a*mfa + b*mfb
315// CastroRadiation is the only code that uses this function.
316template <typename MF>
318FabSetT<MF>::linComb (value_type a, const MF& mfa, int a_comp,
319 value_type b, const MF& mfb, int b_comp,
320 int dcomp, int ncomp, int ngrow)
321{
322 BL_PROFILE("FabSetT<MF>::linComb()");
323 BL_ASSERT(mfa.nGrowVect().allGE(ngrow) && mfb.nGrowVect().allGE(ngrow));
324 BL_ASSERT(mfa.boxArray() == mfb.boxArray());
325 BL_ASSERT(boxArray() != mfa.boxArray());
326
327 MF bdrya(boxArray(),DistributionMap(),ncomp,0,MFInfo());
328 MF bdryb(boxArray(),DistributionMap(),ncomp,0,MFInfo());
329
330 const auto huge = static_cast<value_type>(sizeof(value_type) == 8 ? 1.e200 : 1.e30);
331
332#ifdef AMREX_USE_OMP
333#pragma omp parallel if (Gpu::notInLaunchRegion())
334#endif
335 for (MFIter mfi(bdrya); mfi.isValid(); ++mfi) // tiling is not safe for this BoxArray
336 {
337 const Box& bx = mfi.validbox();
338 auto afab = bdrya.array(mfi);
339 auto bfab = bdryb.array(mfi);
340 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
341 {
342 afab(i,j,k,n) = huge;
343 bfab(i,j,k,n) = huge;
344 });
345 }
346
347 bdrya.ParallelCopy(mfa,a_comp,0,ncomp,ngrow,0);
348 bdryb.ParallelCopy(mfb,b_comp,0,ncomp,ngrow,0);
349
350#ifdef AMREX_USE_OMP
351#pragma omp parallel if (Gpu::notInLaunchRegion())
352#endif
353 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi)
354 {
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);
359 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
360 {
361 dfab(i,j,k,n+dcomp) = a*afab(i,j,k,n) + b*bfab(i,j,k,n);
362 });
363 }
364
365 return *this;
366}
367
368template <typename MF>
369void
370FabSetT<MF>::write (const std::string& name) const
371{
372 if (AsyncOut::UseAsyncOut()) {
373 VisMF::AsyncWrite(m_mf,name);
374 } else {
375 VisMF::Write(m_mf,name);
376 }
377}
378
379template <typename MF>
380void
381FabSetT<MF>::read (const std::string& name)
382{
383 if (m_mf.empty()) {
384 amrex::Abort("FabSetT<MF>::read: not predefined");
385 }
386 VisMF::Read(m_mf,name);
387}
388
389template <typename MF>
390void
392{
395 int ncomp = dst.nComp();
396#ifdef AMREX_USE_OMP
397#pragma omp parallel if (Gpu::notInLaunchRegion())
398#endif
399 for (FabSetIter fsi(dst); fsi.isValid(); ++fsi) {
400 const Box& bx = fsi.validbox();
401 auto const srcfab = src.array(fsi);
402 auto dstfab = dst.array(fsi);
403 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
404 {
405 dstfab(i,j,k,n) = srcfab(i,j,k,n);
406 });
407 }
408}
409
412
413}
414
415#endif /*_FABSET_H_*/
#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