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
17#ifdef AMREX_USE_OMP
18#include <omp.h>
19#endif
20
21#include <limits>
22
23namespace amrex {
24
49template <typename MF>
51{
52 friend class FabSetIter;
53 friend class FluxRegister;
54public:
57
58 //
60 FabSetT () noexcept = default;
61 //
63 FabSetT (const BoxArray& grids, const DistributionMapping& dmap, int ncomp);
64
65 ~FabSetT () = default;
66
67 FabSetT (FabSetT<MF>&& rhs) noexcept = default;
68
69 FabSetT (const FabSetT<MF>& rhs) = delete;
70 FabSetT<MF>& operator= (const FabSetT<MF>& rhs) = delete;
71 FabSetT<MF>& operator= (FabSetT<MF>&& rhs) = delete;
72
73 //
75 void define (const BoxArray& grids, const DistributionMapping& dmap, int ncomp);
76
78 FAB const& operator[] (const MFIter& mfi) const noexcept { return m_mf[mfi]; }
79
81 FAB & operator[] (const MFIter& mfi) noexcept { return m_mf[mfi]; }
82
84 FAB const& operator[] (int i) const noexcept { return m_mf[i]; }
85
87 FAB & operator[] (int i) noexcept { return m_mf[i]; }
88
90 Array4<value_type const> array (const MFIter& mfi) const noexcept { return m_mf.const_array(mfi); }
91
93 Array4<value_type > array (const MFIter& mfi) noexcept { return m_mf.array(mfi); }
94
96 Array4<value_type const> array (int i) const noexcept { return m_mf.const_array(i); }
97
99 Array4<value_type > array (int i) noexcept { return m_mf.array(i); }
100
102 Array4<value_type const> const_array (const MFIter& mfi) const noexcept { return m_mf.const_array(mfi); }
103
105 Array4<value_type const> const_array (int i) const noexcept { return m_mf.const_array(i); }
106
108 MultiArray4<value_type const> arrays () const noexcept { return m_mf.const_arrays(); }
109
111 MultiArray4<value_type > arrays () noexcept { return m_mf.arrays(); }
112
114 MultiArray4<value_type const> const_arrays () const noexcept { return m_mf.const_arrays(); }
115
117 [[nodiscard]] Box fabbox (int K) const noexcept { return m_mf.fabbox(K); }
118
120 [[nodiscard]] int size () const noexcept { return m_mf.size(); }
121
123 [[nodiscard]] const BoxArray& boxArray () const noexcept { return m_mf.boxArray(); }
124
126 [[nodiscard]] const DistributionMapping& DistributionMap () const noexcept
127 { return m_mf.DistributionMap(); }
128
130 [[nodiscard]] MF & multiFab () noexcept { return m_mf; }
131
133 [[nodiscard]] MF const& multiFab () const noexcept { return m_mf; }
134
136 [[nodiscard]] int nComp () const noexcept { return m_mf.nComp(); }
137
139 void clear () { m_mf.clear(); }
140
150 FabSetT<MF>& copyFrom (const FabSetT<MF>& src, int scomp, int dcomp, int ncomp);
151
162 FabSetT<MF>& copyFrom (const MF& src, int ngrow, int scomp, int dcomp, int ncomp,
163 const Periodicity& period = Periodicity::NonPeriodic());
164
173 FabSetT<MF>& plusFrom (const FabSetT<MF>& src, int scomp, int dcomp, int ncomp);
174
186 FabSetT<MF>& plusFrom (const MF& src, int ngrow, int scomp, int dcomp, int ncomp,
187 const Periodicity& period = Periodicity::NonPeriodic(),
188 bool deterministic = false);
189
200 void copyTo (MF& dest, int ngrow, int scomp, int dcomp, int ncomp,
201 const Periodicity& period = Periodicity::NonPeriodic()) const;
202
213 void plusTo (MF& dest, int ngrow, int scomp, int dcomp, int ncomp,
214 const Periodicity& period = Periodicity::NonPeriodic()) const;
215
217 void setVal (value_type val);
218
220 void setVal (value_type val, int comp, int num_comp);
221
233 int scomp, int dcomp, int ncomp);
234
248 FabSetT<MF>& linComb (value_type a, const MF& mfa, int a_comp,
249 value_type b, const MF& mfb, int b_comp,
250 int dcomp, int ncomp, int ngrow);
251
257 void write (const std::string& name) const;
258
264 void read (const std::string& name);
265
272 static void Copy (FabSetT<MF>& dst, const FabSetT<MF>& src);
273
274private:
275 MF m_mf;
276};
277
280 : public MFIter
281{
282public:
283 template <typename MF>
284 explicit FabSetIter (const FabSetT<MF>& fs)
285 : MFIter(fs.m_mf) { }
286};
287
288template <typename MF>
289FabSetT<MF>::FabSetT (const BoxArray& grids, const DistributionMapping& dmap, int ncomp)
290 : m_mf(grids,dmap,ncomp,0)
291{}
292
293template <typename MF>
294void
295FabSetT<MF>::define (const BoxArray& grids, const DistributionMapping& dm, int ncomp)
296{
297 m_mf.define(grids, dm, ncomp, 0);
298}
299
300template <typename MF>
302FabSetT<MF>::copyFrom (const FabSetT<MF>& src, int scomp, int dcomp, int ncomp)
303{
304 if (boxArray() == src.boxArray() && DistributionMap() == src.DistributionMap()) {
305#ifdef AMREX_USE_OMP
306#pragma omp parallel if (Gpu::notInLaunchRegion())
307#endif
308 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi) {
309 const Box& bx = fsi.validbox();
310 auto const srcfab = src.array(fsi);
311 auto dstfab = this->array(fsi);
312 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
313 {
314 dstfab(i,j,k,n+dcomp) = srcfab(i,j,k,n+scomp);
315 });
316 }
317 } else {
318 m_mf.ParallelCopy(src.m_mf,scomp,dcomp,ncomp);
319 }
320 return *this;
321}
322
323template <typename MF>
325FabSetT<MF>::copyFrom (const MF& src, int ngrow, int scomp, int dcomp, int ncomp,
326 const Periodicity& period)
327{
328 BL_ASSERT(boxArray() != src.boxArray());
329 m_mf.ParallelCopy(src,scomp,dcomp,ncomp,ngrow,0,period);
330 return *this;
331}
332
333template <typename MF>
335FabSetT<MF>::plusFrom (const FabSetT<MF>& src, int scomp, int dcomp, int ncomp)
336{
337 if (boxArray() == src.boxArray() && DistributionMap() == src.DistributionMap()) {
338#ifdef AMREX_USE_OMP
339#pragma omp parallel if (Gpu::notInLaunchRegion())
340#endif
341 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi) {
342 const Box& bx = fsi.validbox();
343 auto const srcfab = src.array(fsi);
344 auto dstfab = this->array(fsi);
345 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
346 {
347 dstfab(i,j,k,n+dcomp) += srcfab(i,j,k,n+scomp);
348 });
349 }
350 } else {
351 amrex::Abort("FabSetT<MF>::plusFrom: parallel plusFrom not supported");
352 }
353 return *this;
354}
355
356template <typename MF>
358FabSetT<MF>::plusFrom (const MF& src, int ngrow, int scomp, int dcomp, int ncomp,
359 const Periodicity& period, bool deterministic)
360{
361 BL_ASSERT(boxArray() != src.boxArray());
362 m_mf.ParallelCopy(src,scomp,dcomp,ncomp,IntVect(ngrow),IntVect(0),period,
363 FabArrayBase::ADD,nullptr,deterministic);
364 return *this;
365}
366
367template <typename MF>
368void
369FabSetT<MF>::copyTo (MF& dest, int ngrow, int scomp, int dcomp, int ncomp,
370 const Periodicity& period) const
371{
372 BL_ASSERT(boxArray() != dest.boxArray());
373 dest.ParallelCopy(m_mf,scomp,dcomp,ncomp,0,ngrow,period);
374}
375
376template <typename MF>
377void
378FabSetT<MF>::plusTo (MF& dest, int ngrow, int scomp, int dcomp, int ncomp,
379 const Periodicity& period) const
380{
381 BL_ASSERT(boxArray() != dest.boxArray());
382 dest.ParallelCopy(m_mf,scomp,dcomp,ncomp,0,ngrow,period,FabArrayBase::ADD);
383}
384
385template <typename MF>
386void
388{
389 const int ncomp = nComp();
390#ifdef AMREX_USE_OMP
391#pragma omp parallel if (Gpu::notInLaunchRegion())
392#endif
393 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi) {
394 const Box& bx = fsi.validbox();
395 auto fab = this->array(fsi);
396 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
397 {
398 fab(i,j,k,n) = val;
399 });
400 }
401}
402
403template <typename MF>
404void
405FabSetT<MF>::setVal (value_type val, int comp, int num_comp)
406{
407#ifdef AMREX_USE_OMP
408#pragma omp parallel if (Gpu::notInLaunchRegion())
409#endif
410 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi) {
411 const Box& bx = fsi.validbox();
412 auto fab = this->array(fsi);
413 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, num_comp, i, j, k, n,
414 {
415 fab(i,j,k,n+comp) = val;
416 });
417 }
418}
419
420// Linear combination this := a*this + b*src
421// Note: corresponding fabsets must be commensurate.
422template <typename MF>
425 int scomp, int dcomp, int ncomp)
426{
427 BL_ASSERT(size() == src.size());
428
429#ifdef AMREX_USE_OMP
430#pragma omp parallel if (Gpu::notInLaunchRegion())
431#endif
432 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi)
433 {
434 const Box& bx = fsi.validbox();
435 auto const srcfab = src.array(fsi);
436 auto dstfab = this->array(fsi);
437 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
438 {
439 dstfab(i,j,k,n+dcomp) = a*dstfab(i,j,k,n+dcomp) + b*srcfab(i,j,k,n+scomp);
440 });
441 }
442 return *this;
443}
444
445// Linear combination: this := a*mfa + b*mfb
446// CastroRadiation is the only code that uses this function.
447template <typename MF>
449FabSetT<MF>::linComb (value_type a, const MF& mfa, int a_comp,
450 value_type b, const MF& mfb, int b_comp,
451 int dcomp, int ncomp, int ngrow)
452{
453 BL_PROFILE("FabSetT<MF>::linComb()");
454 BL_ASSERT(mfa.nGrowVect().allGE(ngrow) && mfb.nGrowVect().allGE(ngrow));
455 BL_ASSERT(mfa.boxArray() == mfb.boxArray());
456 BL_ASSERT(boxArray() != mfa.boxArray());
457
458 MF bdrya(boxArray(),DistributionMap(),ncomp,0,MFInfo());
459 MF bdryb(boxArray(),DistributionMap(),ncomp,0,MFInfo());
460
461 const auto huge = static_cast<value_type>(sizeof(value_type) == 8 ? 1.e200 : 1.e30);
462
463#ifdef AMREX_USE_OMP
464#pragma omp parallel if (Gpu::notInLaunchRegion())
465#endif
466 for (MFIter mfi(bdrya); mfi.isValid(); ++mfi) // tiling is not safe for this BoxArray
467 {
468 const Box& bx = mfi.validbox();
469 auto afab = bdrya.array(mfi);
470 auto bfab = bdryb.array(mfi);
471 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
472 {
473 afab(i,j,k,n) = huge;
474 bfab(i,j,k,n) = huge;
475 });
476 }
477
478 bdrya.ParallelCopy(mfa,a_comp,0,ncomp,ngrow,0);
479 bdryb.ParallelCopy(mfb,b_comp,0,ncomp,ngrow,0);
480
481#ifdef AMREX_USE_OMP
482#pragma omp parallel if (Gpu::notInLaunchRegion())
483#endif
484 for (FabSetIter fsi(*this); fsi.isValid(); ++fsi)
485 {
486 const Box& bx = fsi.validbox();
487 auto const afab = bdrya.array(fsi);
488 auto const bfab = bdryb.array(fsi);
489 auto dfab = this->array(fsi);
490 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
491 {
492 dfab(i,j,k,n+dcomp) = a*afab(i,j,k,n) + b*bfab(i,j,k,n);
493 });
494 }
495
496 return *this;
497}
498
499template <typename MF>
500void
501FabSetT<MF>::write (const std::string& name) const
502{
503 if (AsyncOut::UseAsyncOut()) {
504 VisMF::AsyncWrite(m_mf,name);
505 } else {
506 VisMF::Write(m_mf,name);
507 }
508}
509
510template <typename MF>
511void
512FabSetT<MF>::read (const std::string& name)
513{
514 if (m_mf.empty()) {
515 amrex::Abort("FabSetT<MF>::read: not predefined");
516 }
517 VisMF::Read(m_mf,name);
518}
519
520template <typename MF>
521void
523{
526 int ncomp = dst.nComp();
527#ifdef AMREX_USE_OMP
528#pragma omp parallel if (Gpu::notInLaunchRegion())
529#endif
530 for (FabSetIter fsi(dst); fsi.isValid(); ++fsi) {
531 const Box& bx = fsi.validbox();
532 auto const srcfab = src.array(fsi);
533 auto dstfab = dst.array(fsi);
534 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
535 {
536 dstfab(i,j,k,n) = srcfab(i,j,k,n);
537 });
538 }
539}
540
543
544}
545
546#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:568
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
@ ADD
Definition AMReX_FabArrayBase.H:394
Convenience iterator that reuses MFIter semantics for FabSet traversal.
Definition AMReX_FabSet.H:281
FabSetIter(const FabSetT< MF > &fs)
Definition AMReX_FabSet.H:284
A FabSet is a group of FArrayBox's. The grouping is designed specifically to represent regions along ...
Definition AMReX_FabSet.H:51
Array4< value_type const > const_array(const MFIter &mfi) const noexcept
Return a const Array4 view for iterator mfi.
Definition AMReX_FabSet.H:102
MultiArray4< value_type const > arrays() const noexcept
Return multi-fab const Array4 views to every FAB.
Definition AMReX_FabSet.H:108
int nComp() const noexcept
Number of components carried by each FAB.
Definition AMReX_FabSet.H:136
Array4< value_type const > const_array(int i) const noexcept
Return a const Array4 view for FAB index i.
Definition AMReX_FabSet.H:105
MF & multiFab() noexcept
Mutable access to the underlying MultiFab.
Definition AMReX_FabSet.H:130
Array4< value_type const > array(int i) const noexcept
Return an Array4 view for FAB index i (const).
Definition AMReX_FabSet.H:96
static void Copy(FabSetT< MF > &dst, const FabSetT< MF > &src)
Local copy function.
Definition AMReX_FabSet.H:522
FabSetT< MF > & copyFrom(const FabSetT< MF > &src, int scomp, int dcomp, int ncomp)
Copy components from another FabSet with identical layout.
Definition AMReX_FabSet.H:302
FabSetT< MF > & plusFrom(const FabSetT< MF > &src, int scomp, int dcomp, int ncomp)
Accumulate data from another FabSet with identical layout.
Definition AMReX_FabSet.H:335
void read(const std::string &name)
Read (used for reading from checkpoint).
Definition AMReX_FabSet.H:512
void clear()
Release storage and reset the internal MultiFab.
Definition AMReX_FabSet.H:139
FAB const & operator[](const MFIter &mfi) const noexcept
Access the FAB referenced by iterator mfi (const).
Definition AMReX_FabSet.H:78
void plusTo(MF &dest, int ngrow, int scomp, int dcomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic()) const
Add this boundary data back into a MultiFab.
Definition AMReX_FabSet.H:378
MultiArray4< value_type > arrays() noexcept
Return multi-fab mutable Array4 views to every FAB.
Definition AMReX_FabSet.H:111
FabSetT< MF > & linComb(value_type a, value_type b, const FabSetT< MF > &src, int scomp, int dcomp, int ncomp)
Linear combination: this := a*this + b*src.
Definition AMReX_FabSet.H:424
Box fabbox(int K) const noexcept
Get the bounding Box of FAB index K.
Definition AMReX_FabSet.H:117
Array4< value_type > array(int i) noexcept
Return an Array4 view for FAB index i (mutable).
Definition AMReX_FabSet.H:99
MultiArray4< value_type const > const_arrays() const noexcept
Return multi-fab const Array4 views (alias of arrays()).
Definition AMReX_FabSet.H:114
const BoxArray & boxArray() const noexcept
Return the BoxArray defining FAB locations.
Definition AMReX_FabSet.H:123
void write(const std::string &name) const
Write (used for writing to checkpoint).
Definition AMReX_FabSet.H:501
FabSetT() noexcept=default
The default constructor – you must later call define().
typename FabDataType< MF >::value_type value_type
Definition AMReX_FabSet.H:55
void setVal(value_type val)
Fill every component of every FAB with scalar value val.
Definition AMReX_FabSet.H:387
int size() const noexcept
Return number of FABs stored in this set.
Definition AMReX_FabSet.H:120
Array4< value_type const > array(const MFIter &mfi) const noexcept
Return an Array4 view for iterator mfi (const).
Definition AMReX_FabSet.H:90
typename FabDataType< MF >::fab_type FAB
Definition AMReX_FabSet.H:56
Array4< value_type > array(const MFIter &mfi) noexcept
Return an Array4 view for iterator mfi (mutable).
Definition AMReX_FabSet.H:93
const DistributionMapping & DistributionMap() const noexcept
Return the DistributionMapping describing FAB ownership.
Definition AMReX_FabSet.H:126
void define(const BoxArray &grids, const DistributionMapping &dmap, int ncomp)
Define a FabSetT<MF> constructed via default constructor on grids/dmap with ncomp components.
Definition AMReX_FabSet.H:295
MF const & multiFab() const noexcept
Const access to the underlying MultiFab.
Definition AMReX_FabSet.H:133
void copyTo(MF &dest, int ngrow, int scomp, int dcomp, int ncomp, const Periodicity &period=Periodicity::NonPeriodic()) const
Copy data from this boundary set back into a MultiFab.
Definition AMReX_FabSet.H:369
Flux Register.
Definition AMReX_FluxRegister.H:25
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
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:2851
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2866
bool match(const BoxArray &x, const BoxArray &y)
Note that two BoxArrays that match are not necessarily equal.
Definition AMReX_BoxArray.cpp:1934
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2861
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Definition AMReX_FabDataType.H:9
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_FabArray.H:153