Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_AuxBoundaryData.H
Go to the documentation of this file.
1
2#ifndef AMREX_AuxBoundaryData_H_
3#define AMREX_AuxBoundaryData_H_
4#include <AMReX_Config.H>
5
6#include <AMReX_Geometry.H>
7#include <AMReX_MultiFab.H>
8
9namespace amrex {
11{
12public:
13
14 AuxBoundaryData () noexcept = default;
15
16 AuxBoundaryData (const BoxArray& ba,
17 int n_grow,
18 int n_comp,
19 const Geometry& geom);
20
21 ~AuxBoundaryData () = default;
22
24 AuxBoundaryData& operator= (AuxBoundaryData&& rhs) = default;
25
27 AuxBoundaryData& operator= (const AuxBoundaryData& rhs) = delete;
28
29 void copyTo (MultiFab& destmf,
30 int src_comp,
31 int dst_comp,
32 int num_comp) const;
33
34 void copyFrom (const MultiFab& srcmf,
35 int src_comp,
36 int dst_comp,
37 int num_comp,
38 int src_ng = 0);
39
40 size_t size () const noexcept
41 {
42 BL_ASSERT(!m_empty); BL_ASSERT(m_initialized); return m_fabs.size();
43 }
44
45 void copy (const AuxBoundaryData& src,
46 int src_comp,
47 int dst_comp,
48 int num_comp);
49
50 void initialize (const BoxArray& ba,
51 int n_grow,
52 int n_comp,
53 const Geometry& geom);
54
55 const BoxArray& equivBoxArray () const noexcept
56 {
57 BL_ASSERT(!m_empty); BL_ASSERT(m_initialized); return m_fabs.boxArray();
58 }
59
60 void setVal (Real r) { BL_ASSERT(m_initialized); if (!m_empty) { m_fabs.setVal(r); } }
61
62 const DistributionMapping& DistributionMap () const noexcept
63 {
64 BL_ASSERT(!m_empty); BL_ASSERT(m_initialized); return m_fabs.DistributionMap();
65 }
66
67 FArrayBox& operator[] (const MFIter& mfi) noexcept
68 {
69 BL_ASSERT(!m_empty); BL_ASSERT(m_initialized); return m_fabs[mfi];
70 }
71 const FArrayBox& operator[] (const MFIter& mfi) const noexcept
72 {
73 BL_ASSERT(!m_empty); BL_ASSERT(m_initialized); return m_fabs[mfi];
74 }
75
76 int nGrow () const noexcept { BL_ASSERT(m_initialized); return m_ngrow; }
77
78 int nComp () const noexcept
79 {
80 BL_ASSERT(!m_empty); BL_ASSERT(m_initialized); return m_fabs.nComp();
81 }
82
83 bool isEmpty () const noexcept { return m_empty; }
84
85private:
86
87 MultiFab m_fabs;
88 int m_ngrow{0};
89 bool m_empty{false};
90 bool m_initialized{false};
91};
92
93}
94
95#endif
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
Definition AMReX_AuxBoundaryData.H:11
int nComp() const noexcept
Definition AMReX_AuxBoundaryData.H:78
const BoxArray & equivBoxArray() const noexcept
Definition AMReX_AuxBoundaryData.H:55
bool isEmpty() const noexcept
Definition AMReX_AuxBoundaryData.H:83
size_t size() const noexcept
Definition AMReX_AuxBoundaryData.H:40
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_AuxBoundaryData.H:62
AuxBoundaryData() noexcept=default
int nGrow() const noexcept
Definition AMReX_AuxBoundaryData.H:76
void copyFrom(const MultiFab &srcmf, int src_comp, int dst_comp, int num_comp, int src_ng=0)
Definition AMReX_AuxBoundaryData.cpp:132
FArrayBox & operator[](const MFIter &mfi) noexcept
Definition AMReX_AuxBoundaryData.H:67
void copyTo(MultiFab &destmf, int src_comp, int dst_comp, int num_comp) const
Definition AMReX_AuxBoundaryData.cpp:118
void initialize(const BoxArray &ba, int n_grow, int n_comp, const Geometry &geom)
Definition AMReX_AuxBoundaryData.cpp:47
void setVal(Real r)
Definition AMReX_AuxBoundaryData.H:60
void copy(const AuxBoundaryData &src, int src_comp, int dst_comp, int num_comp)
Definition AMReX_AuxBoundaryData.cpp:22
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
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
int size() const noexcept
Return the number of FABs in the FabArray.
Definition AMReX_FabArrayBase.H:110
const DistributionMapping & DistributionMap() const noexcept
Return constant reference to associated DistributionMapping.
Definition AMReX_FabArrayBase.H:131
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition AMReX_FabArrayBase.H:83
const BoxArray & boxArray() const noexcept
Return a constant reference to the BoxArray that defines the valid region associated with this FabArr...
Definition AMReX_FabArrayBase.H:95
void setVal(value_type val)
Set all components in the entire region of each FAB to val.
Definition AMReX_FabArray.H:2652
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:63
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
Definition AMReX_Amr.cpp:49